Skip to content

Commit

Permalink
Merge pull request #794 from aritorto/missedCommentPR793
Browse files Browse the repository at this point in the history
[refactor] LevelCartesianIndexMapper for PolyhedralGrid to eliminate repeated code
  • Loading branch information
lisajulia authored Dec 20, 2024
2 parents b7ef2ef + 17a2306 commit e693927
Show file tree
Hide file tree
Showing 3 changed files with 17 additions and 35 deletions.
2 changes: 1 addition & 1 deletion opm/grid/cpgrid/LevelCartesianIndexMapper.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -52,7 +52,7 @@ template<>
class LevelCartesianIndexMapper<Dune::CpGrid>
{
public:
static const int dimension = 3 ;
static constexpr int dimension = 3 ;

explicit LevelCartesianIndexMapper(const Dune::CpGrid& grid) : grid_{ &grid }
{}
Expand Down
47 changes: 15 additions & 32 deletions opm/grid/polyhedralgrid/levelcartesianindexmapper.hh
Original file line number Diff line number Diff line change
Expand Up @@ -50,72 +50,55 @@ namespace Opm
// Interface class to access the local Cartesian grid of each level grid (when refinement).
// Further documentation in opm/grid/common/LevelCartesianIndexMapper.hpp
//
// Adapter Design Pattern: In this case, LevelCartesianIndexMapper uses the Object Adapter variant, where it holds an instance
// (here, a std::unique_ptr) of CartesianIndexMapper, the wrapped type. The goal is to provide a standardized interface, allowing
// incompatible functionality (such as Cartesian indexing in the context of refinement that may not be supported - yet -for all
// grid types, like CpGrid) to integrate smoothly within the existing conventions.
//
// Specialization for PolyhedralGrid
template<int dim, int dimworld, typename coord_t>
class LevelCartesianIndexMapper<Dune::PolyhedralGrid< dim, dimworld, coord_t >>
{
typedef Dune::PolyhedralGrid< dim, dimworld, coord_t > Grid;
using Grid = Dune::PolyhedralGrid< dim, dimworld, coord_t >;
public:
static const int dimension = 3 ;
static constexpr int dimension = 3 ;

explicit LevelCartesianIndexMapper(const Grid& grid) : grid_{ &grid }
explicit LevelCartesianIndexMapper(const Dune::CartesianIndexMapper<Grid>& cartesian_index_mapper)
: cartesianIndexMapper_{std::make_unique<Dune::CartesianIndexMapper<Grid>>(cartesian_index_mapper)}
{}

const std::array<int,3>& cartesianDimensions(int level) const
{
throwIfLevelPositive(level);
return grid_->logicalCartesianSize();
return cartesianIndexMapper_->logicalCartesianSize();
}

int cartesianSize(int level) const
{
throwIfLevelPositive(level);
return computeCartesianSize(0);
return cartesianIndexMapper_->cartesianSize();
}

int compressedSize(int level) const
{
throwIfLevelPositive(level);
return grid_->size(0);
return cartesianIndexMapper_->compressedSize();
}

int cartesianIndex( const int compressedElementIndex, const int level) const
{
throwIfLevelPositive(level);
assert( compressedElementIndex >= 0 && compressedElementIndex < compressedSize(0) );
return grid_->globalCell()[ compressedElementIndex ];
return cartesianIndexMapper_->cartesianIndex(compressedElementIndex);
}

void cartesianCoordinate(const int compressedElementIndex, std::array<int,dimension>& coords, int level) const
{
throwIfLevelPositive(level);

int gc = cartesianIndex( compressedElementIndex, 0);
auto cartesianDimensions = grid_->logicalCartesianSize();
if( dimension >=2 )
{
for( int d=0; d<dimension-2; ++d )
{
coords[d] = gc % cartesianDimensions[d]; gc /= cartesianDimensions[d];
}

coords[dimension-2] = gc % cartesianDimensions[dimension-2];
coords[dimension-1] = gc / cartesianDimensions[dimension-1];
}
else
coords[ 0 ] = gc ;
cartesianIndexMapper_->cartesianCoordinate(compressedElementIndex, coords);
}

private:
const Grid* grid_;

int computeCartesianSize(int level) const
{
int size = cartesianDimensions(level)[ 0 ];
for( int d=1; d<dimension; ++d )
size *= cartesianDimensions(level)[ d ];
return size;
}
std::unique_ptr<Dune::CartesianIndexMapper<Grid>> cartesianIndexMapper_;

void throwIfLevelPositive(int level) const
{
Expand Down
3 changes: 1 addition & 2 deletions tests/test_lookupdata_polyhedral.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -91,15 +91,14 @@ void lookup_check(const Dune::PolyhedralGrid<3,3>& grid)

std::vector<double> fake_feature_double(grid.size(0), 0.);
std::iota(fake_feature_double.begin(), fake_feature_double.end(), .5);

const Opm::LevelCartesianIndexMapper< Dune::PolyhedralGrid<3,3> > levelCartMapp(grid);

const auto leaf_view = grid.leafGridView();
using GridView = std::remove_cv_t< typename std::remove_reference<decltype(grid.leafGridView())>::type>;
// LookUpData
const Opm::LookUpData<Dune::PolyhedralGrid<3,3>, GridView> lookUpData(leaf_view, false);
// LookUpCartesianData
const Dune::CartesianIndexMapper<Dune::PolyhedralGrid<3,3>> cartMapper(grid);
const Opm::LevelCartesianIndexMapper< Dune::PolyhedralGrid<3,3> > levelCartMapp(cartMapper);
const Opm::LookUpCartesianData<Dune::PolyhedralGrid<3,3>, GridView> lookUpCartesianData(leaf_view, cartMapper, false);
// Mapper
const Dune::MultipleCodimMultipleGeomTypeMapper<GridView> mapper(grid.leafGridView(), Dune::mcmgElementLayout());
Expand Down

0 comments on commit e693927

Please sign in to comment.