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
32#include <opm/models/common/multiphasebaseproperties.hh>
33
36
37#include <array>
38#include <functional>
39#include <string>
40#include <tuple>
41#include <unordered_set>
42
43namespace Opm {
44template <class TypeTag>
46}
47
48namespace Opm::Properties {
49
50namespace TTag {
52 using InheritsFrom = std::tuple<FlowBaseVanguard>;
53};
54}
55
56// declare the properties
57template<class TypeTag>
58struct Vanguard<TypeTag, TTag::PolyhedralGridVanguard> {
60};
61template<class TypeTag>
62struct Grid<TypeTag, TTag::PolyhedralGridVanguard> {
63 using type = Dune::PolyhedralGrid<3, 3>;
64};
65template<class TypeTag>
66struct EquilGrid<TypeTag, TTag::PolyhedralGridVanguard> {
67 using type = GetPropType<TypeTag, Properties::Grid>;
68};
69
70} // namespace Opm::Properties
71
72namespace Opm {
73
81template <class TypeTag>
83{
84 friend class FlowBaseVanguard<TypeTag>;
85 using ParentType = FlowBaseVanguard<TypeTag>;
86
87 using ElementMapper = GetPropType<TypeTag, Properties::ElementMapper>;
88 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
89 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
90
91public:
92 using Grid = GetPropType<TypeTag, Properties::Grid>;
93 using EquilGrid = GetPropType<TypeTag, Properties::EquilGrid>;
94 using GridView = GetPropType<TypeTag, Properties::GridView>;
97 static constexpr int dimension = Grid::dimension;
98 static constexpr int dimensionworld = Grid::dimensionworld;
99
100private:
101 using GridPointer = Grid*;
102 using EquilGridPointer = EquilGrid*;
103
104public:
106 CartesianIndexMapper, Scalar>;
107
108 PolyhedralGridVanguard(Simulator& simulator)
109 : FlowBaseVanguard<TypeTag>(simulator)
110 , simulator_(simulator)
111 {
113 // add a copy in standard vector format to fullfill new interface
114 const int* globalcellorg = this->grid().globalCell();
115 int num_cells = this->gridView().size(0);
116 globalcell_.resize(num_cells);
117 for(int i=0; i < num_cells; ++i){
118 globalcell_[i] = globalcellorg[i];
119 }
120 }
121
126 { return *grid_; }
127
131 const Grid& grid() const
132 { return *grid_; }
133
143 const EquilGrid& equilGrid() const
144 { return *grid_; }
145
154 { /* do nothing: The EQUIL grid is the simulation grid! */ }
155
162 { /* do nothing: PolyhedralGrid is not parallel! */
163 }
164
170 { return *cartesianIndexMapper_; }
171
179 { return *cartesianIndexMapper_; }
180
181 const std::vector<int>& globalCell()
182 {
183 return globalcell_;
184 }
185
186 unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const {
187 return elemIndex;
188 }
189
190 unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const {
191 return elemIndex;
192 }
193
200 {
201 }
202
203 std::unordered_set<std::string> defunctWellNames() const
204 { return defunctWellNames_; }
205
207 {
208 return simulator_.problem().eclTransmissibilities();
209 }
210
218 std::function<std::array<double,FlowBaseVanguard<TypeTag>::dimensionworld>(int)>
220 {
221 return this->cellCentroids_(this->cartesianIndexMapper(), false);
222 }
223
224 std::vector<int> cellPartition() const
225 {
226 // not required for this type of grid yet (only from bdaBridge??)
227 return {};
228 }
229protected:
231 {
232 grid_ = std::make_unique<Grid>(this->eclState().getInputGrid(), this->eclState().fieldProps().porv(true));
233 cartesianIndexMapper_ = std::make_unique<CartesianIndexMapper>(*grid_);
234 this->updateGridView_();
236 this->updateCellDepths_();
237 }
238
240 {
241 // not handling the removal of completions for this type of grid yet.
242 }
243
244 Simulator& simulator_;
245
246 std::unique_ptr<Grid> grid_;
247 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
248 //CartesianIndexMapperPointer cartesianIndexMapper_;
249
250 std::unordered_set<std::string> defunctWellNames_;
251 std::vector<int> globalcell_;
252};
253
254} // namespace Opm
255
256#endif // OPM_POLYHEDRAL_GRID_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
Helper class for grid instantiation of ECL file-format using problems.
Definition: FlowBaseVanguard.hpp:216
void updateCartesianToCompressedMapping_()
Definition: FlowBaseVanguard.hpp:526
void updateCellDepths_()
Definition: FlowBaseVanguard.hpp:548
void callImplementationInit()
Definition: FlowBaseVanguard.hpp:516
GetPropType< TypeTag, Properties::GridView > GridView
Definition: FlowBaseVanguard.hpp:227
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:496
const EclipseState & eclState() const
Return a reference to the internalized ECL deck.
Definition: FlowGenericVanguard.hpp:120
Helper class for grid instantiation of ECL file-format using problems.
Definition: PolyhedralGridVanguard.hpp:83
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
Definition: PolyhedralGridVanguard.hpp:153
std::unique_ptr< Grid > grid_
Definition: PolyhedralGridVanguard.hpp:246
unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const
Definition: PolyhedralGridVanguard.hpp:186
const CartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
Definition: PolyhedralGridVanguard.hpp:178
void createGrids_()
Definition: PolyhedralGridVanguard.hpp:230
static constexpr int dimensionworld
Definition: PolyhedralGridVanguard.hpp:98
std::unordered_set< std::string > defunctWellNames_
Definition: PolyhedralGridVanguard.hpp:250
void releaseGlobalTransmissibilities()
Free the memory occupied by the global transmissibility object.
Definition: PolyhedralGridVanguard.hpp:199
const std::vector< int > & globalCell()
Definition: PolyhedralGridVanguard.hpp:181
void loadBalance()
Distribute the simulation grid over multiple processes.
Definition: PolyhedralGridVanguard.hpp:161
Simulator & simulator_
Definition: PolyhedralGridVanguard.hpp:244
Dune::CartesianIndexMapper< Grid > CartesianIndexMapper
Definition: PolyhedralGridVanguard.hpp:95
static constexpr int dimension
Definition: PolyhedralGridVanguard.hpp:97
const TransmissibilityType & globalTransmissibility() const
Definition: PolyhedralGridVanguard.hpp:206
std::vector< int > cellPartition() const
Definition: PolyhedralGridVanguard.hpp:224
const EquilGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
Definition: PolyhedralGridVanguard.hpp:143
GetPropType< TypeTag, Properties::GridView > GridView
Definition: PolyhedralGridVanguard.hpp:94
GetPropType< TypeTag, Properties::EquilGrid > EquilGrid
Definition: PolyhedralGridVanguard.hpp:93
std::vector< int > globalcell_
Definition: PolyhedralGridVanguard.hpp:251
PolyhedralGridVanguard(Simulator &simulator)
Definition: PolyhedralGridVanguard.hpp:108
std::unordered_set< std::string > defunctWellNames() const
Definition: PolyhedralGridVanguard.hpp:203
void filterConnections_()
Definition: PolyhedralGridVanguard.hpp:239
Grid & grid()
Return a reference to the simulation grid.
Definition: PolyhedralGridVanguard.hpp:125
GetPropType< TypeTag, Properties::Grid > Grid
Definition: PolyhedralGridVanguard.hpp:92
std::unique_ptr< CartesianIndexMapper > cartesianIndexMapper_
Definition: PolyhedralGridVanguard.hpp:247
const Grid & grid() const
Return a reference to the simulation grid.
Definition: PolyhedralGridVanguard.hpp:131
unsigned int gridIdxToEquilGridIdx(unsigned int elemIndex) const
Definition: PolyhedralGridVanguard.hpp:190
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
Definition: PolyhedralGridVanguard.hpp:169
std::function< std::array< double, FlowBaseVanguard< TypeTag >::dimensionworld >(int)> cellCentroids() const
Get function to query cell centroids for a distributed grid.
Definition: PolyhedralGridVanguard.hpp:219
Definition: Transmissibility.hpp:54
Definition: AluGridVanguard.hpp:57
Definition: BlackoilPhases.hpp:27
GetPropType< TypeTag, Properties::Grid > type
Definition: PolyhedralGridVanguard.hpp:67
Definition: FlowBaseVanguard.hpp:64
Dune::PolyhedralGrid< 3, 3 > type
Definition: PolyhedralGridVanguard.hpp:63
Definition: PolyhedralGridVanguard.hpp:51
std::tuple< FlowBaseVanguard > InheritsFrom
Definition: PolyhedralGridVanguard.hpp:52