Split disjoint columns in extractColumns(). Add test case.

This commit is contained in:
Atgeirr Flø Rasmussen 2012-05-11 11:20:18 +02:00
parent 891696e333
commit 5b2e99ec4c
2 changed files with 165 additions and 36 deletions

View File

@ -5,29 +5,48 @@
namespace Opm {
/// Helper struct for extractColumn
/// Compares the underlying k-index
struct ExtractColumnCompare
{
ExtractColumnCompare(const UnstructuredGrid& g)
: grid(g)
{
// empty
}
namespace {
bool operator()(const int i, const int j)
{
// Extract k-index
int index_i = grid.global_cell ? grid.global_cell[i] : i;
int k_i = index_i / grid.cartdims[0] / grid.cartdims[1];
int index_j = grid.global_cell ? grid.global_cell[j] : j;
int k_j = index_j / grid.cartdims[0] / grid.cartdims[1];
/// Helper struct for extractColumn
/// Compares the underlying k-index
struct ExtractColumnCompare
{
ExtractColumnCompare(const UnstructuredGrid& g)
: grid(g)
{
// empty
}
return k_i < k_j;
}
bool operator()(const int i, const int j)
{
// Extract k-index
int index_i = grid.global_cell ? grid.global_cell[i] : i;
int k_i = index_i / grid.cartdims[0] / grid.cartdims[1];
int index_j = grid.global_cell ? grid.global_cell[j] : j;
int k_j = index_j / grid.cartdims[0] / grid.cartdims[1];
return k_i < k_j;
}
const UnstructuredGrid& grid;
};
/// Neighbourhood query.
/// \return true if two cells are neighbours.
bool neighbours(const UnstructuredGrid& grid, const int c0, const int c1)
{
for (int hf = grid.cell_facepos[c0]; hf < grid.cell_facepos[c0 + 1]; ++hf) {
const int f = grid.cell_faces[hf];
if (grid.face_cells[2*f] == c1 || grid.face_cells[2*f+1] == c1) {
return true;
}
}
return false;
}
} // anonymous namespace
const UnstructuredGrid& grid;
};
/// Extract each column of the grid.
/// \note Assumes the pillars of the grid are all vertically aligned.
@ -64,6 +83,40 @@ inline void extractColumn( const UnstructuredGrid& grid, std::vector<std::vector
for (int col = 0; col < num_cols; ++col) {
std::sort(columns[col].begin(), columns[col].end(), ExtractColumnCompare(grid));
}
// At this point, a column may contain multiple disjoint sets of cells.
// We must split these columns into connected parts.
std::vector< std::vector<int> > new_columns;
for (int col = 0; col < num_cols; ++col) {
const int colsz = columns[col].size();
int first_of_col = 0;
for (int k = 1; k < colsz; ++k) {
const int c0 = columns[col][k-1];
const int c1 = columns[col][k];
if (!neighbours(grid, c0, c1)) {
// Must split. Move the cells [first_of_col, ... , k-1] to
// a new column, known to be connected.
new_columns.push_back(std::vector<int>());
new_columns.back().assign(columns[col].begin() + first_of_col, columns[col].begin() + k);
// The working column now starts with index k.
first_of_col = k;
}
}
if (first_of_col != 0) {
// The column was split, the working part should be
// the entire column. We erase the cells before first_of_col.
// (Could be more efficient if we instead chop off end.)
columns[col].erase(columns[col].begin(), columns[col].begin() + first_of_col);
}
}
// Must tack on the new columns to complete the set.
const int num_cols_all = num_cols + new_columns.size();
columns.resize(num_cols_all);
for (int col = num_cols; col < num_cols_all; ++col) {
columns[col].swap(new_columns[col - num_cols]);
}
}
} // namespace Opm

View File

@ -8,6 +8,7 @@
#include <boost/test/unit_test.hpp>
#include <opm/core/ColumnExtract.hpp>
#include <opm/core/GridManager.hpp>
#include <opm/core/eclipse/EclipseGridParser.hpp>
BOOST_AUTO_TEST_CASE(SingleColumnTest)
{
@ -16,21 +17,14 @@ BOOST_AUTO_TEST_CASE(SingleColumnTest)
const int size_x = 1, size_y = 1, size_z = 10;
GridManager manager(size_x, size_y, size_z);
// We do our own numbering
UnstructuredGrid* grid = const_cast<UnstructuredGrid*>(manager.c_grid());
grid->global_cell = (int*)malloc(sizeof(int) * size_x * size_y * size_z);
for(int i = 0; i < size_x * size_y * size_z; ++i) {
grid->global_cell[i] = i;
}
std::vector<std::vector<int> > columns;
extractColumn(*manager.c_grid(), columns);
std::vector<int> correct_answer;
for(int i = 0; i < 10; i++) {
for (int i = 0; i < 10; ++i) {
correct_answer.push_back(i);
}
BOOST_CHECK_EQUAL_COLLECTIONS(correct_answer.begin(), correct_answer.end(),
columns[0].begin(), columns[0].end());
@ -43,13 +37,6 @@ BOOST_AUTO_TEST_CASE(FourByFourColumnTest)
using namespace Opm;
GridManager manager(size_x, size_y, size_z);
// We do our own numbering
UnstructuredGrid* grid = const_cast<UnstructuredGrid*>(manager.c_grid());
grid->global_cell = (int*)malloc(sizeof(int) * size_x * size_y * size_z);
for(int i = 0; i < size_x * size_y * size_z; ++i) {
grid->global_cell[i] = i;
}
std::vector<std::vector<int> > columns;
extractColumn(*manager.c_grid(), columns);
@ -67,3 +54,92 @@ BOOST_AUTO_TEST_CASE(FourByFourColumnTest)
}
}
BOOST_AUTO_TEST_CASE(DisjointColumn)
{
std::string grdecl =
"SPECGRID \n"
"3 3 3 1 F \n"
"/ \n"
" \n"
"COORD \n"
"0.0 0 0 0.0 0 3.0 \n"
"1.0 0 0 1.0 0 3.0 \n"
"2.0 0 0 2.0 0 3.0 \n"
"3.0 0 0 3.0 0 3.0 \n"
" \n"
"0.0 1.0 0 0.0 1.0 3.0 \n"
"1.0 1.0 0 1.0 1.0 3.0 \n"
"2.0 1.0 0 2.0 1.0 3.0 \n"
"3.0 1.0 0 3.0 1.0 3.0 \n"
" \n"
"0.0 2.0 0 0.0 2.0 3.0 \n"
"1.0 2.0 0 1.0 2.0 3.0 \n"
"2.0 2.0 0 2.0 2.0 3.0 \n"
"3.0 2.0 0 3.0 2.0 3.0 \n"
" \n"
"0.0 3.0 0 0.0 3.0 3.0 \n"
"1.0 3.0 0 1.0 3.0 3.0 \n"
"2.0 3.0 0 2.0 3.0 3.0 \n"
"3.0 3.0 0 3.0 3.0 3.0 \n"
"/ \n"
" \n"
"ZCORN \n"
"36*0.0 \n"
"72*1.0 \n"
"72*2.0 \n"
"36*3.0 \n"
"/ \n"
" \n"
"ACTNUM \n"
"13*1 \n"
"0 \n"
"13*1 \n"
"/ \n"
"\n";
using namespace Opm;
std::istringstream is(grdecl);
EclipseGridParser deck;
deck.read(is);
GridManager manager(deck);
std::vector<std::vector<int> > columns;
extractColumn(*manager.c_grid(), columns);
// for (size_t i = 0; i < columns.size(); ++i) {
// for (size_t j = 0; j < columns[i].size(); ++j) {
// std::cout << columns[i][j] << ' ';
// }
// std::cout << '\n';
// }
// std::cout << std::endl;
const int correct[10][3] = { { 0, 1, 2 },
{ 3, 4, 5 },
{ 6, 7, 8 },
{ 9, 10, 11 },
{ 13, -1, -1 },
{ 14, 15, 16 },
{ 17, 18, 19 },
{ 20, 21, 22 },
{ 23, 24, 25 },
{ 12, -1, -1 } };
std::vector<std::vector<int> > correct_answer;
correct_answer.resize(10);
for(int i = 0; i < 10; i++) {
for(int j = 0; j < 3; j++) {
correct_answer[i].push_back(correct[i][j]);
}
}
correct_answer[4].resize(1);
correct_answer[9].resize(1);
BOOST_CHECK_EQUAL(columns.size(), 10);
for(int i = 0; i < 10; i++) {
BOOST_CHECK_EQUAL_COLLECTIONS(correct_answer[i].begin(), correct_answer[i].end(),
columns[i].begin(), columns[i].end());
}
}