PolyhedralGridVanguard.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
27#ifndef OPM_POLYHEDRAL_GRID_VANGUARD_HPP
28#define OPM_POLYHEDRAL_GRID_VANGUARD_HPP
29
30#include <opm/grid/polyhedralgrid.hh>
31#include <opm/grid/polyhedralgrid/levelcartesianindexmapper.hh>
32
34
37
38#include <array>
39#include <functional>
40#include <string>
41#include <tuple>
42#include <unordered_set>
43
44namespace Opm {
45template <class TypeTag>
47}
48
49namespace Opm::Properties {
50
51namespace TTag {
53 using InheritsFrom = std::tuple<FlowBaseVanguard>;
54};
55}
56
57// declare the properties
58template<class TypeTag>
59struct Vanguard<TypeTag, TTag::PolyhedralGridVanguard> {
61};
62template<class TypeTag>
63struct Grid<TypeTag, TTag::PolyhedralGridVanguard> {
64 using type = Dune::PolyhedralGrid<3, 3>;
65};
66template<class TypeTag>
67struct EquilGrid<TypeTag, TTag::PolyhedralGridVanguard> {
69};
70
71} // namespace Opm::Properties
72
73namespace Opm {
74
82template <class TypeTag>
84{
85 friend class FlowBaseVanguard<TypeTag>;
87
91
92public:
99 static constexpr int dimension = Grid::dimension;
100 static constexpr int dimensionworld = Grid::dimensionworld;
101
102private:
103 using GridPointer = Grid*;
104 using EquilGridPointer = EquilGrid*;
105
106public:
108 CartesianIndexMapper, Scalar>;
109
110 explicit PolyhedralGridVanguard(Simulator& simulator)
111 : FlowBaseVanguard<TypeTag>(simulator)
112 , simulator_(simulator)
113 {
115 // add a copy in standard vector format to fullfill new interface
116 const int* globalcellorg = this->grid().globalCell();
117 int num_cells = this->gridView().size(0);
118 globalcell_.resize(num_cells);
119 for(int i=0; i < num_cells; ++i){
120 globalcell_[i] = globalcellorg[i];
121 }
122 }
123
128 { return *grid_; }
129
133 const Grid& grid() const
134 { return *grid_; }
135
145 const EquilGrid& equilGrid() const
146 { return *grid_; }
147
156 { /* do nothing: The EQUIL grid is the simulation grid! */ }
157
164 { /* do nothing: PolyhedralGrid is not parallel! */
165 }
166
167 void addLgrs()
168 { /* do nothing: PolyhedralGrid with LGRs not supported yet! */
169 }
170
176 { return *cartesianIndexMapper_; }
177
185
193 { return *cartesianIndexMapper_; }
194
195 const std::vector<int>& globalCell()
196 {
197 return globalcell_;
198 }
199
200 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const {
201 return elemIndex;
202 }
203
204 unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const {
205 return elemIndex;
206 }
207
214 {
215 }
216
217 std::unordered_set<std::string> defunctWellNames() const
218 { return defunctWellNames_; }
219
221 {
222 return simulator_.problem().eclTransmissibilities();
223 }
224
232 std::function<std::array<double,FlowBaseVanguard<TypeTag>::dimensionworld>(int)>
234 {
235 return this->cellCentroids_(this->cartesianIndexMapper(), false);
236 }
237
238 std::vector<int> cellPartition() const
239 {
240 // not required for this type of grid yet (only from bdaBridge??)
241 return {};
242 }
243
244protected:
246 {
247 this->grid_ = std::make_unique<Grid>
248 (this->eclState().getInputGrid(),
249 this->eclState().fieldProps().porv(true));
250
252 std::make_unique<CartesianIndexMapper>(*this->grid_);
253
254 this->updateGridView_();
256 this->updateCellDepths_();
257 }
258
260 {
261 // not handling the removal of completions for this type of grid yet.
262 }
263
264 Simulator& simulator_;
265
266 std::unique_ptr<Grid> grid_;
267 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
268 //CartesianIndexMapperPointer cartesianIndexMapper_;
269
270 std::unordered_set<std::string> defunctWellNames_;
271 std::vector<int> globalcell_;
272};
273
274} // namespace Opm
275
276#endif // OPM_POLYHEDRAL_GRID_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
Provides the base class for most (all?) simulator vanguards.
Definition: basevanguard.hh:50
void updateGridView_()
Definition: basevanguard.hh:132
const GridView & gridView() const
Returns a reference to the grid view to be used.
Definition: basevanguard.hh:70
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:84
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:333
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:355
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:321
std::function< std::array< double, dimensionworld >(int)> cellCentroids_(const CartMapper &cartMapper, const bool &isCpGrid) const
Get function to query cell centroids for a distributed grid.
Definition: FlowBaseVanguard.hpp:301
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:167
Definition: RelpermDiagnostics.hpp:31
Helper class for grid instantiation of ECL file-format using problems.
Definition: PolyhedralGridVanguard.hpp:84
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition: PolyhedralGridVanguard.hpp:155
std::unique_ptr< Grid > grid_
Definition: PolyhedralGridVanguard.hpp:266
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
Definition: PolyhedralGridVanguard.hpp:200
const CartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition: PolyhedralGridVanguard.hpp:192
void createGrids_()
Definition: PolyhedralGridVanguard.hpp:245
static constexpr int dimensionworld
Definition: PolyhedralGridVanguard.hpp:100
std::unordered_set< std::string > defunctWellNames_
Definition: PolyhedralGridVanguard.hpp:270
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition: PolyhedralGridVanguard.hpp:213
const std::vector< int > & globalCell()
Definition: PolyhedralGridVanguard.hpp:195
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition: PolyhedralGridVanguard.hpp:163
Simulator & simulator_
Definition: PolyhedralGridVanguard.hpp:264
Dune::CartesianIndexMapper< Grid > CartesianIndexMapper
Definition: PolyhedralGridVanguard.hpp:96
static constexpr int dimension
Definition: PolyhedralGridVanguard.hpp:99
const TransmissibilityType & globalTransmissibility() const
Definition: PolyhedralGridVanguard.hpp:220
std::vector< int > cellPartition() const
Definition: PolyhedralGridVanguard.hpp:238
const EquilGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition: PolyhedralGridVanguard.hpp:145
GetPropType< TypeTag, Properties::GridView > GridView
Definition: PolyhedralGridVanguard.hpp:95
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: PolyhedralGridVanguard.hpp:94
std::vector< int > globalcell_
Definition: PolyhedralGridVanguard.hpp:271
PolyhedralGridVanguard(Simulator &simulator)
Definition: PolyhedralGridVanguard.hpp:110
std::unordered_set< std::string > defunctWellNames() const
Definition: PolyhedralGridVanguard.hpp:217
void filterConnections_()
Definition: PolyhedralGridVanguard.hpp:259
Grid & grid()
Return a reference to the simulation grid.
Definition: PolyhedralGridVanguard.hpp:127
Opm::LevelCartesianIndexMapper< Grid > LevelCartesianIndexMapper
Definition: PolyhedralGridVanguard.hpp:97
GetPropType< TypeTag, Properties::Grid > Grid
Definition: PolyhedralGridVanguard.hpp:93
std::unique_ptr< CartesianIndexMapper > cartesianIndexMapper_
Definition: PolyhedralGridVanguard.hpp:267
const Grid & grid() const
Return a reference to the simulation grid.
Definition: PolyhedralGridVanguard.hpp:133
void addLgrs()
Definition: PolyhedralGridVanguard.hpp:167
unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const
Definition: PolyhedralGridVanguard.hpp:204
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: PolyhedralGridVanguard.hpp:175
const LevelCartesianIndexMapper levelCartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: PolyhedralGridVanguard.hpp:183
std::function< std::array< double, FlowBaseVanguard< TypeTag >::dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition: PolyhedralGridVanguard.hpp:233
Definition: Transmissibility.hpp:54
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilmodel.hh:79
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
GetPropType< TypeTag, Properties::Grid > type
Definition: PolyhedralGridVanguard.hpp:68
Definition: FlowBaseVanguard.hpp:70
Dune::PolyhedralGrid< 3, 3 > type
Definition: PolyhedralGridVanguard.hpp:64
The type of the DUNE grid.
Definition: basicproperties.hh:100
Definition: PolyhedralGridVanguard.hpp:52
std::tuple< FlowBaseVanguard > InheritsFrom
Definition: PolyhedralGridVanguard.hpp:53
Property which provides a Vanguard (manages grids)
Definition: basicproperties.hh:96