Add orderCounterClockwise() and test.

This commit is contained in:
Atgeirr Flø Rasmussen
2014-11-21 08:54:26 +01:00
parent 9bf36c37d4
commit b26bfc51f6
3 changed files with 76 additions and 2 deletions

View File

@@ -20,6 +20,8 @@
#include <opm/core/grid/GridUtilities.hpp>
#include <set>
#include <vector>
#include <cmath>
#include <algorithm>
namespace Opm
{
@@ -72,4 +74,51 @@ namespace Opm
// 3. Done. Return.
return cell_nb;
}
/// For each cell, order the (cell) neighbours counterclockwise.
/// \param[in] grid A 2d grid object.
/// \param[in, out] nb A cell-cell neighbourhood table, such as from vertexNeighbours().
void orderCounterClockwise(const UnstructuredGrid& grid,
SparseTable<int>& nb)
{
if (grid.dimensions != 2) {
OPM_THROW(std::logic_error, "Cannot use orderCounterClockwise in " << grid.dimensions << " dimensions.");
}
const int num_cells = grid.number_of_cells;
if (nb.size() != num_cells) {
OPM_THROW(std::logic_error, "Inconsistent arguments for orderCounterClockwise().");
}
// For each cell, compute each neighbour's angle with the x axis,
// sort that to find the correct permutation of the neighbours.
typedef std::pair<double, int> AngleAndPos;
std::vector<AngleAndPos> angle_and_pos;
std::vector<int> original;
for (int cell = 0; cell < num_cells; ++cell) {
const int num_nb = nb[cell].size();
angle_and_pos.clear();
angle_and_pos.resize(num_nb);
for (int ii = 0; ii < num_nb; ++ii) {
const int cell2 = nb[cell][ii];
const double v[2] = { grid.cell_centroids[2*cell2] - grid.cell_centroids[2*cell],
grid.cell_centroids[2*cell2 + 1] - grid.cell_centroids[2*cell + 1] };
double angle = std::acos(v[0]/std::sqrt(v[0]*v[0] + v[1]*v[1]));
if (v[1] < 0.0) {
angle = 2*M_PI - angle;
}
angle_and_pos[ii] = std::make_pair(angle, ii);
}
original.assign(nb[cell].begin(), nb[cell].end());
std::sort(angle_and_pos.begin(), angle_and_pos.end());
for (int ii = 0; ii < num_nb; ++ii) {
nb[cell][ii] = original[angle_and_pos[ii].second];
}
}
}
} // namespace Opm

View File

@@ -31,6 +31,12 @@ namespace Opm
/// \return A table of neighbour cell-indices by cell.
SparseTable<int> vertexNeighbours(const UnstructuredGrid& grid);
/// For each cell, order the (cell) neighbours counterclockwise.
/// \param[in] grid A 2d grid object.
/// \param[in, out] nb A cell-cell neighbourhood table, such as from vertexNeighbours().
void orderCounterClockwise(const UnstructuredGrid& grid,
SparseTable<int>& nb);
} // namespace Opm
#endif // OPM_GRIDUTILITIES_HEADER_INCLUDED

View File

@@ -33,7 +33,7 @@
using namespace Opm;
BOOST_AUTO_TEST_CASE(cartesian_2d)
BOOST_AUTO_TEST_CASE(cartesian_2d_vertexNeighbours)
{
const GridManager gm(2, 2);
const UnstructuredGrid& grid = *gm.c_grid();
@@ -47,7 +47,7 @@ BOOST_AUTO_TEST_CASE(cartesian_2d)
BOOST_CHECK(vnb == truth);
}
BOOST_AUTO_TEST_CASE(cartesian_3d)
BOOST_AUTO_TEST_CASE(cartesian_3d_vertexNeighbours)
{
const GridManager gm(3, 2, 2);
const UnstructuredGrid& grid = *gm.c_grid();
@@ -60,3 +60,22 @@ BOOST_AUTO_TEST_CASE(cartesian_3d)
const int nb[n] = { 1, 3, 4, 6, 7, 9, 10 };
BOOST_CHECK_EQUAL_COLLECTIONS(vnb[0].begin(), vnb[0].end(), nb, nb + n);
}
BOOST_AUTO_TEST_CASE(cartesian_2d_orderCounterClockwise)
{
const GridManager gm(2, 2);
const UnstructuredGrid& grid = *gm.c_grid();
SparseTable<int> vnb = vertexNeighbours(grid);
orderCounterClockwise(grid, vnb);
BOOST_REQUIRE(!vnb.empty());
const int num_elem = 12;
const int elem[num_elem] = { 1, 3, 2, 3, 2, 0, 3, 0, 1, 2, 0, 1 };
const int num_rows = 4;
const int rowsizes[num_rows] = { 3, 3, 3, 3 };
const SparseTable<int> truth(elem, elem + num_elem, rowsizes, rowsizes + num_rows);
BOOST_CHECK(vnb == truth);
for (int c = 0; c < num_rows; ++c) {
BOOST_CHECK_EQUAL_COLLECTIONS(vnb[c].begin(), vnb[c].end(), truth[c].begin(), truth[c].end());
}
}