opm-grid
Loading...
Searching...
No Matches
LevelCartesianIndexMapper.hpp
1//===========================================================================
2//
3// File: LevelCartesianIndexMapper.hpp
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_CPGRIDLEVELCARTESIANINDEXMAPPER_HH
35#define OPM_CPGRIDLEVELCARTESIANINDEXMAPPER_HH
36
37#include <opm/grid/common/LevelCartesianIndexMapper.hpp>
38#include <opm/grid/CpGrid.hpp>
39
40namespace Dune
41{
42class CpGrid;
43}
44
45namespace Opm
46{
47// Interface class to access the local Cartesian grid of each level grid (when refinement).
48// Further documentation in opm/grid/common/LevelCartesianIndexMapper.hpp
49//
50// Specialization for CpGrid
51template<>
52class LevelCartesianIndexMapper<Dune::CpGrid>
53{
54public:
55 static constexpr int dimension = 3 ;
56
57 explicit LevelCartesianIndexMapper(const Dune::CpGrid& grid) : grid_{ &grid }
58 {}
59
60 LevelCartesianIndexMapper() = delete;
61
62 const std::array<int,3>& cartesianDimensions(int level) const
63 {
64 validLevel(level);
65 return grid_->currentData()[level]->logicalCartesianSize();
66 }
67
68 int cartesianSize(int level) const
69 {
70 validLevel(level);
71 return computeCartesianSize(level);
72 }
73
74 int compressedSize(int level) const
75 {
76 validLevel(level);
77 return grid_->currentData()[level]->size(0);
78 }
79
80 int cartesianIndex( const int compressedElementIndex, const int level) const
81 {
82 validLevel(level);
83 assert( compressedElementIndex >= 0 && compressedElementIndex < grid_->currentData()[level]->size(0) );
84 return grid_->currentData()[level]->globalCell()[compressedElementIndex];
85 }
86
87 void cartesianCoordinate(const int compressedElementIndexOnLevel, std::array<int,dimension>& coordsOnLevel, int level) const
88 {
89 validLevel(level);
90 grid_->currentData()[level]->getIJK( compressedElementIndexOnLevel, coordsOnLevel);
91 }
92
93private:
94 const Dune::CpGrid* grid_;
95
96 int computeCartesianSize(int level) const
97 {
98 int size = cartesianDimensions(level)[ 0 ];
99 for( int d=1; d<dimension; ++d )
100 size *= cartesianDimensions(level)[ d ];
101 return size;
102 }
103
104 void validLevel(int level) const
105 {
106 if ((level < 0) || (level > grid_->maxLevel())) {
107 throw std::invalid_argument("Invalid level.\n");
108 }
109 }
110};
111
112}
113
114#endif
[ provides Dune::Grid ]
Definition CpGrid.hpp:203
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68