Correct column extract
This commit is contained in:
parent
565430c2e9
commit
6a0df4a752
@ -1,25 +1,55 @@
|
||||
#include <opm/core/grid.h>
|
||||
#include <vector>
|
||||
#include <map>
|
||||
#include <algorithm>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
/// Extract each column of the grid.
|
||||
/// \note Assumes the pillars of the grid are all vertically aligned.
|
||||
/// \param grid The grid from which to extract the columns.
|
||||
/// \param columns will for each i + cartgrim[0]*j contain the k values contained in the column
|
||||
/// centered at (i, j).
|
||||
void extractColumn( const UnstructuredGrid& grid, std::vector<std::vector<int> >& columns )
|
||||
/// Helper struct for extractColumn
|
||||
/// Compares the underlying k-index
|
||||
struct ExtractColumnCompare
|
||||
{
|
||||
ExtractColumnCompare(const UnstructuredGrid& grid)
|
||||
: grid(grid)
|
||||
{
|
||||
const int* dims = grid.cartdims;
|
||||
columns.resize(dims[0]*dims[1]);
|
||||
for (int i = 0; i < dims[0]; ++i) {
|
||||
for (int j = 0; j < dims[1]; ++j) {
|
||||
int plane_index = i + j*dims[0];
|
||||
for (int k = 0; k < dims[2]; ++k) {
|
||||
columns[plane_index].push_back(plane_index + k*dims[0]*dims[1]);
|
||||
}
|
||||
}
|
||||
}
|
||||
// empty
|
||||
}
|
||||
|
||||
bool operator()(const int i, const int j)
|
||||
{
|
||||
// Extract k-index
|
||||
int k_i = grid.global_cell[i] / grid.cartdims[0] / grid.cartdims[1];
|
||||
int k_j = grid.global_cell[j] / grid.cartdims[0] / grid.cartdims[1];
|
||||
|
||||
return k_i <= k_j;
|
||||
}
|
||||
|
||||
const UnstructuredGrid& grid;
|
||||
};
|
||||
|
||||
/// Extract each column of the grid.
|
||||
/// \note Assumes the pillars of the grid are all vertically aligned.
|
||||
/// \param grid The grid from which to extract the columns.
|
||||
/// \param columns will for each i + cartgrim[0]*j contain the k values contained in the column
|
||||
/// centered at (i, j).
|
||||
void extractColumn( const UnstructuredGrid& grid, std::map<int, std::vector<int> >& columns )
|
||||
{
|
||||
const int* dims = grid.cartdims;
|
||||
const int* global = grid.global_cell;
|
||||
for (int i = 0; i < grid.number_of_cells; ++i) {
|
||||
// Extract Cartesian coordinates
|
||||
int index = grid.global_cell[i];
|
||||
int i_cart = index % dims[0];
|
||||
int k_cart = index / dims[0] / dims[1];
|
||||
int j_cart = (index - k_cart*dims[0]*dims[1])/ dims[0];
|
||||
|
||||
// Finally add it to the map
|
||||
columns[i_cart + j_cart*dims[0]].push_back(i);
|
||||
}
|
||||
|
||||
for (std::map<int, std::vector<int> >::iterator it = columns.begin(); it != columns.end(); ++it) {
|
||||
std::sort(it->second.begin(), it->second.end(), ExtractColumnCompare(grid));
|
||||
}
|
||||
}
|
||||
|
||||
} // namespace Opm
|
||||
|
@ -11,11 +11,19 @@
|
||||
|
||||
BOOST_AUTO_TEST_CASE(SingleColumnTest)
|
||||
{
|
||||
const int size_x = 1, size_y = 1, size_z = 10;
|
||||
using namespace Opm;
|
||||
|
||||
const int size_x = 1, size_y = 1, size_z = 10;
|
||||
GridManager manager(size_x, size_y, size_z);
|
||||
|
||||
std::vector<std::vector<int> > columns;
|
||||
// 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::map<int, std::vector<int> > columns;
|
||||
extractColumn(*manager.c_grid(), columns);
|
||||
|
||||
std::vector<int> correct_answer;
|
||||
@ -34,7 +42,14 @@ BOOST_AUTO_TEST_CASE(FourByFourColumnTest)
|
||||
using namespace Opm;
|
||||
GridManager manager(size_x, size_y, size_z);
|
||||
|
||||
std::vector<std::vector<int> > columns;
|
||||
// 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::map<int, std::vector<int> > columns;
|
||||
extractColumn(*manager.c_grid(), columns);
|
||||
|
||||
std::vector<std::vector<int> > correct_answer;
|
||||
|
Loading…
Reference in New Issue
Block a user