GenericCpGridVanguard.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_GENERIC_CPGRID_VANGUARD_HPP
28#define OPM_GENERIC_CPGRID_VANGUARD_HPP
29
30#include <opm/grid/CpGrid.hpp>
31
33
34#include <functional>
35#include <memory>
36#include <optional>
37#include <vector>
38
39#if HAVE_MPI
40#include <filesystem>
41#endif
42
43namespace Opm {
44 class EclipseState;
45 class Schedule;
46 class Well;
47 class ParallelEclipseState;
48 class LgrCollection;
49}
50
51namespace Opm {
52
53#if HAVE_MPI
54namespace details {
56{
57public:
58 explicit MPIPartitionFromFile(const std::filesystem::path& partitionFile)
59 : partitionFile_(partitionFile)
60 {}
61
62 std::vector<int> operator()(const Dune::CpGrid& grid) const;
63
64private:
65 std::filesystem::path partitionFile_{};
66};
67
68} // namespace Opm::details
69#endif // HAVE_MPI
70
71
75extern std::optional<std::function<std::vector<int> (const Dune::CpGrid&)>> externalLoadBalancer;
76
77template<class ElementMapper, class GridView, class Scalar>
79protected:
81 using Element = typename GridView::template Codim<0>::Entity;
82
83public:
85
86 virtual ~GenericCpGridVanguard() = default;
87
91 Dune::CpGrid& grid()
92 { return *grid_; }
93
97 const Dune::CpGrid& grid() const
98 { return *grid_; }
99
109 const Dune::CpGrid& equilGrid() const;
110
119
123 static void setExternalLoadBalancer(const std::function<std::vector<int> (const Dune::CpGrid&)>& loadBalancer)
124 {
125 externalLoadBalancer = loadBalancer;
126 }
127
133
138
139 const std::vector<int>& cellPartition() const
140 {
141 return this->cell_part_;
142 }
143
144protected:
150#if HAVE_MPI
151 void doLoadBalance_(const Dune::EdgeWeightMethod edgeWeightsMethod,
152 const bool ownersFirst,
153 const bool serialPartitioning,
154 const bool enableDistributedWells,
155 const double zoltanImbalanceTol,
156 const GridView& gridView,
157 const Schedule& schedule,
158 EclipseState& eclState,
160 const int numJacobiBlocks);
161
162 void distributeFieldProps_(EclipseState& eclState);
163
164private:
165 std::vector<double> extractFaceTrans(const GridView& gridView) const;
166
167 void distributeGrid(const Dune::EdgeWeightMethod edgeWeightsMethod,
168 const bool ownersFirst,
169 const bool serialPartitioning,
170 const bool enableDistributedWells,
171 const double zoltanImbalanceTol,
172 const bool loadBalancerSet,
173 const std::vector<double>& faceTrans,
174 const std::vector<Well>& wells,
175 EclipseState& eclState,
177
178 void distributeGrid(const Dune::EdgeWeightMethod edgeWeightsMethod,
179 const bool ownersFirst,
180 const bool serialPartitioning,
181 const bool enableDistributedWells,
182 const double zoltanImbalanceTol,
183 const bool loadBalancerSet,
184 const std::vector<double>& faceTrans,
185 const std::vector<Well>& wells,
186 ParallelEclipseState* eclState,
188
189protected:
190 virtual const std::string& zoltanParams() const = 0;
191
192#endif // HAVE_MPI
193
195
196 void doCreateGrids_(EclipseState& eclState);
197 void addLgrsUpdateLeafView(const LgrCollection& lgrCollection, const int lgrsSize);
198
199 virtual void allocTrans() = 0;
200 virtual double getTransmissibility(unsigned I, unsigned J) const = 0;
201
202 // removing some connection located in inactive grid cells
203 void doFilterConnections_(Schedule& schedule);
204
205 Scalar computeCellThickness(const Element& element) const;
206
207 std::unique_ptr<Dune::CpGrid> grid_;
208 std::unique_ptr<Dune::CpGrid> equilGrid_;
209 std::unique_ptr<CartesianIndexMapper> cartesianIndexMapper_;
210 std::unique_ptr<CartesianIndexMapper> equilCartesianIndexMapper_;
211
213 std::vector<int> cell_part_{};
214};
215
216} // namespace Opm
217
218#endif // OPM_GENERIC_CPGRID_VANGUARD_HPP
Definition: CollectDataOnIORank.hpp:49
std::vector< std::pair< std::string, bool > > ParallelWellStruct
Definition: FlowGenericVanguard.hpp:62
Definition: GenericCpGridVanguard.hpp:78
const CartesianIndexMapper & cartesianIndexMapper() const
Returns the object which maps a global element index of the simulation grid to the corresponding elem...
void doLoadBalance_(const Dune::EdgeWeightMethod edgeWeightsMethod, const bool ownersFirst, const bool serialPartitioning, const bool enableDistributedWells, const double zoltanImbalanceTol, const GridView &gridView, const Schedule &schedule, EclipseState &eclState, FlowGenericVanguard::ParallelWellStruct &parallelWells, const int numJacobiBlocks)
Distribute the simulation grid over multiple processes.
void releaseEquilGrid()
Indicates that the initial condition has been computed and the memory used by the EQUIL grid can be r...
const Dune::CpGrid & equilGrid() const
Returns a refefence to the grid which should be used by the EQUIL initialization code.
std::unique_ptr< Dune::CpGrid > equilGrid_
Definition: GenericCpGridVanguard.hpp:208
virtual double getTransmissibility(unsigned I, unsigned J) const =0
void doCreateGrids_(EclipseState &eclState)
const CartesianIndexMapper & equilCartesianIndexMapper() const
Returns mapper from compressed to cartesian indices for the EQUIL grid.
std::vector< int > cell_part_
Definition: GenericCpGridVanguard.hpp:213
virtual void allocTrans()=0
int mpiRank
Definition: GenericCpGridVanguard.hpp:212
const std::vector< int > & cellPartition() const
Definition: GenericCpGridVanguard.hpp:139
const Dune::CpGrid & grid() const
Return a reference to the simulation grid.
Definition: GenericCpGridVanguard.hpp:97
void distributeFieldProps_(EclipseState &eclState)
std::unique_ptr< CartesianIndexMapper > equilCartesianIndexMapper_
Definition: GenericCpGridVanguard.hpp:210
static void setExternalLoadBalancer(const std::function< std::vector< int >(const Dune::CpGrid &)> &loadBalancer)
Sets a function that returns external load balancing information when passed the grid.
Definition: GenericCpGridVanguard.hpp:123
Dune::CpGrid & grid()
Return a reference to the simulation grid.
Definition: GenericCpGridVanguard.hpp:91
typename GridView::template Codim< 0 >::Entity Element
Definition: GenericCpGridVanguard.hpp:81
virtual const std::string & zoltanParams() const =0
std::unique_ptr< Dune::CpGrid > grid_
Definition: GenericCpGridVanguard.hpp:207
void addLgrsUpdateLeafView(const LgrCollection &lgrCollection, const int lgrsSize)
std::unique_ptr< CartesianIndexMapper > cartesianIndexMapper_
Definition: GenericCpGridVanguard.hpp:209
Scalar computeCellThickness(const Element &element) const
void doFilterConnections_(Schedule &schedule)
virtual ~GenericCpGridVanguard()=default
Parallel frontend to the EclipseState.
Definition: ParallelEclipseState.hpp:147
Definition: GenericCpGridVanguard.hpp:56
MPIPartitionFromFile(const std::filesystem::path &partitionFile)
Definition: GenericCpGridVanguard.hpp:58
std::vector< int > operator()(const Dune::CpGrid &grid) const
Definition: BlackoilPhases.hpp:27
std::optional< std::function< std::vector< int >(const Dune::CpGrid &)> > externalLoadBalancer
optional functor returning external load balancing information