levelcartesianindexmapper.hh
Go to the documentation of this file.
1//===========================================================================
2//
3// File: levelcartesianindexmapper.hh
4//
5// Created: Tue October 01 09:44:00 2024
6//
7// Author(s): Antonella Ritorto <antonella.ritorto@opm-op.com>
8//
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2024 Equinor ASA.
18
19 This file is part of The Open Porous Media project (OPM).
20
21 OPM is free software: you can redistribute it and/or modify
22 it under the terms of the GNU General Public License as published by
23 the Free Software Foundation, either version 3 of the License, or
24 (at your option) any later version.
25
26 OPM is distributed in the hope that it will be useful,
27 but WITHOUT ANY WARRANTY; without even the implied warranty of
28 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
29 GNU General Public License for more details.
30
31 You should have received a copy of the GNU General Public License
32 along with OPM. If not, see <http://www.gnu.org/licenses/>.
33*/
34#ifndef OPM_POLYHEDRALGRIDLEVELCARTESIANINDEXMAPPER_HH
35#define OPM_POLYHEDRALGRIDLEVELCARTESIANINDEXMAPPER_HH
36
39
40
41namespace Dune
42{
43template<int dim, int dimworld, typename coord_t>
44class PolyhedralGrid;
45
46}
47
48namespace Opm
49{
50// Interface class to access the local Cartesian grid of each level grid (when refinement).
51// Further documentation in opm/grid/common/LevelCartesianIndexMapper.hpp
52//
53// Adapter Design Pattern: In this case, LevelCartesianIndexMapper uses the Object Adapter variant, where it holds an instance
54// (here, a std::unique_ptr) of CartesianIndexMapper, the wrapped type. The goal is to provide a standardized interface, allowing
55// incompatible functionality (such as Cartesian indexing in the context of refinement that may not be supported - yet -for all
56// grid types, like CpGrid) to integrate smoothly within the existing conventions.
57//
58// Specialization for PolyhedralGrid
59template<int dim, int dimworld, typename coord_t>
60class LevelCartesianIndexMapper<Dune::PolyhedralGrid< dim, dimworld, coord_t >>
61{
63public:
64 static constexpr int dimension = 3 ;
65
66 explicit LevelCartesianIndexMapper(const Dune::CartesianIndexMapper<Grid>& cartesian_index_mapper)
67 : cartesianIndexMapper_{std::make_unique<Dune::CartesianIndexMapper<Grid>>(cartesian_index_mapper)}
68 {}
69
70 const std::array<int,3>& cartesianDimensions(int level) const
71 {
72 throwIfLevelPositive(level);
73 return cartesianIndexMapper_->logicalCartesianSize();
74 }
75
76 int cartesianSize(int level) const
77 {
78 throwIfLevelPositive(level);
79 return cartesianIndexMapper_->cartesianSize();
80 }
81
82 int compressedSize(int level) const
83 {
84 throwIfLevelPositive(level);
85 return cartesianIndexMapper_->compressedSize();
86 }
87
88 int cartesianIndex( const int compressedElementIndex, const int level) const
89 {
90 throwIfLevelPositive(level);
91 return cartesianIndexMapper_->cartesianIndex(compressedElementIndex);
92 }
93
94 void cartesianCoordinate(const int compressedElementIndex, std::array<int,dimension>& coords, int level) const
95 {
96 throwIfLevelPositive(level);
97 cartesianIndexMapper_->cartesianCoordinate(compressedElementIndex, coords);
98 }
99
100private:
101 std::unique_ptr<Dune::CartesianIndexMapper<Grid>> cartesianIndexMapper_;
102
103 void throwIfLevelPositive(int level) const
104 {
105 if (level) {
106 throw std::invalid_argument("Invalid level.\n");
107 }
108 }
109};
110
111}
112
113#endif
Interface class to access the logical Cartesian grid as used in industry standard simulator decks.
Definition: common/CartesianIndexMapper.hpp:16
identical grid wrapper
Definition: grid.hh:159
const std::array< int, 3 > & cartesianDimensions(int level) const
Definition: levelcartesianindexmapper.hh:70
void cartesianCoordinate(const int compressedElementIndex, std::array< int, dimension > &coords, int level) const
Definition: levelcartesianindexmapper.hh:94
LevelCartesianIndexMapper(const Dune::CartesianIndexMapper< Grid > &cartesian_index_mapper)
Definition: levelcartesianindexmapper.hh:66
int cartesianSize(int level) const
Definition: levelcartesianindexmapper.hh:76
int cartesianIndex(const int compressedElementIndex, const int level) const
Definition: levelcartesianindexmapper.hh:88
int compressedSize(int level) const
Definition: levelcartesianindexmapper.hh:82
Definition: common/LevelCartesianIndexMapper.hpp:45
static constexpr int dimension
Definition: common/LevelCartesianIndexMapper.hpp:48
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26
STL namespace.