GenericTemperatureModel_impl.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*/
28#ifndef OPM_GENERIC_TEMPERATURE_MODEL_IMPL_HPP
29#define OPM_GENERIC_TEMPERATURE_MODEL_IMPL_HPP
30
31#include <dune/istl/operators.hh>
32#include <dune/istl/solvers.hh>
33#include <dune/istl/schwarz.hh>
34#include <dune/istl/preconditioners.hh>
35#include <dune/istl/schwarz.hh>
36
37#include <opm/common/OpmLog/OpmLog.hpp>
38
39#include <opm/grid/CpGrid.hpp>
40
41#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
42#include <opm/input/eclipse/Schedule/Well/Well.hpp>
43
45
50
51#include <fmt/format.h>
52
53#include <array>
54#include <functional>
55#include <memory>
56#include <set>
57#include <stdexcept>
58#include <string>
59
60namespace Opm {
61
62#if HAVE_MPI
63template<class M, class V>
65{
66 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
67 using EnergyOperator = Dune::OverlappingSchwarzOperator<M, V, V, Comm>;
69};
70
71template<class Vector, class Grid, class Matrix>
72std::tuple<std::unique_ptr<Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
73 Dune::OwnerOverlapCopyCommunication<int,int>>>,
74 std::unique_ptr<typename EnergySolverSelector<Matrix,Vector>::type>>
75createParallelFlexibleSolver(const Grid&, const Matrix&, const PropertyTree&)
76{
77 OPM_THROW(std::logic_error, "Grid not supported for parallel Temperatures.");
78 return {nullptr, nullptr};
79}
81template<class Vector, class Matrix>
82std::tuple<std::unique_ptr<Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
83 Dune::OwnerOverlapCopyCommunication<int,int>>>,
84 std::unique_ptr<typename EnergySolverSelector<Matrix,Vector>::type>>
85createParallelFlexibleSolver(const Dune::CpGrid& grid, const Matrix& M, const PropertyTree& prm)
86{
87 using EnergyOperator = Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
88 Dune::OwnerOverlapCopyCommunication<int,int>>;
89 using EnergySolver = Dune::FlexibleSolver<EnergyOperator>;
90 const auto& cellComm = grid.cellCommunication();
91 auto op = std::make_unique<EnergyOperator>(M, cellComm);
92 auto dummyWeights = [](){ return Vector();};
93 return {std::move(op), std::make_unique<EnergySolver>(*op, cellComm, prm, dummyWeights, 0)};
94}
95#endif
96
97template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
99GenericTemperatureModel(const GridView& gridView,
100 const EclipseState& eclState,
101 const CartesianIndexMapper& cartMapper,
102 const DofMapper& dofMapper)
103 : gridView_(gridView)
104 , eclState_(eclState)
105 , cartMapper_(cartMapper)
106 , dofMapper_(dofMapper)
107{
108}
109
110
111template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
113doInit(std::size_t numGridDof)
114{
115 doTemp_ = eclState_.getSimulationConfig().isTemp();
116
117 temperature_.resize(numGridDof);
118 energyVector_.resize(numGridDof);
119 // allocate matrix for storing the Jacobian of the temperature residual
120 energyMatrix_ = std::make_unique<EnergyMatrix>(numGridDof, numGridDof, EnergyMatrix::random);
121
122 // find the sparsity pattern of the temperature matrix
123 using NeighborSet = std::set<unsigned>;
124 std::vector<NeighborSet> neighbors(numGridDof);
125
126 Stencil stencil(gridView_, dofMapper_);
127 for (const auto& elem : elements(gridView_)) {
128 stencil.update(elem);
129
130 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
131 unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
132
133 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
134 unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
135 neighbors[myIdx].insert(neighborIdx);
136 }
137 }
138 }
139
140 // allocate space for the rows of the matrix
141 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
142 energyMatrix_->setrowsize(dofIdx, neighbors[dofIdx].size());
143 }
144 energyMatrix_->endrowsizes();
145
146 // fill the rows with indices. each degree of freedom talks to
147 // all of its neighbors. (it also talks to itself since
148 // degrees of freedom are sometimes quite egocentric.)
149 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
150 typename NeighborSet::iterator nIt = neighbors[dofIdx].begin();
151 typename NeighborSet::iterator nEndIt = neighbors[dofIdx].end();
152 for (; nIt != nEndIt; ++nIt) {
153 energyMatrix_->addindex(dofIdx, *nIt);
154 }
155 }
156 energyMatrix_->endindices();
157
158 maxTempChange_ = Parameters::Get<Parameters::MaxTemperatureChange<Scalar>>();
159}
160
161template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
164{
165 x = 0.0;
166 Scalar tolerance = 1e-2;
167 int maxIter = 100;
168
169 int verbosity = 0;
170 PropertyTree prm;
171 prm.put("maxiter", maxIter);
172 prm.put("tol", tolerance);
173 prm.put("verbosity", verbosity);
174 prm.put("solver", std::string("bicgstab"));
175 prm.put("preconditioner.type", std::string("ParOverILU0"));
176
177#if HAVE_MPI
178 if(gridView_.grid().comm().size() > 1)
179 {
180 auto [energyOperator, solver] =
181 createParallelFlexibleSolver<EnergyVector>(gridView_.grid(), M, prm);
182 (void) energyOperator;
183
185 solver->apply(x, b, result);
186
187 // return the result of the solver
188 return result.converged;
189 }
190 else
191 {
192#endif
193 using EnergySolver = Dune::BiCGSTABSolver<EnergyVector>;
194 using EnergyOperator = Dune::MatrixAdapter<EnergyMatrix,EnergyVector,EnergyVector>;
195 using EnergyScalarProduct = Dune::SeqScalarProduct<EnergyVector>;
196 using EnergyPreconditioner = Dune::SeqILU< EnergyMatrix,EnergyVector,EnergyVector>;
197
198 EnergyOperator energyOperator(M);
199 EnergyScalarProduct energyScalarProduct;
200 EnergyPreconditioner energyPreconditioner(M, 0, 1); // results in ILU0
201
202 EnergySolver solver (energyOperator, energyScalarProduct,
203 energyPreconditioner, tolerance, maxIter,
204 verbosity);
205
207 solver.apply(x, b, result);
208
209 // return the result of the solver
210 return result.converged;
211#if HAVE_MPI
212 }
213#endif
214}
215
216} // namespace Opm
217
218#endif // OPM_GENERIC_TEMPERATURE_MODEL_IMPL_HPP
Definition: CollectDataOnIORank.hpp:49
Definition: FlexibleSolver.hpp:45
GenericTemperatureModel(const GridView &gridView, const EclipseState &eclState, const CartesianIndexMapper &cartMapper, const DofMapper &dofMapper)
Definition: GenericTemperatureModel_impl.hpp:99
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > EnergyMatrix
Definition: GenericTemperatureModel.hpp:57
bool linearSolve_(const EnergyMatrix &M, EnergyVector &x, EnergyVector &b)
Definition: GenericTemperatureModel_impl.hpp:163
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > EnergyVector
Definition: GenericTemperatureModel.hpp:58
void doInit(std::size_t numGridDof)
Initialize all internal data structures needed by the temperature module.
Definition: GenericTemperatureModel_impl.hpp:113
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
void put(const std::string &key, const T &data)
Definition: blackoilbioeffectsmodules.hh:43
Dune::InverseOperatorResult InverseOperatorResult
Definition: GpuBridge.hpp:32
std::tuple< std::unique_ptr< Dune::OverlappingSchwarzOperator< Matrix, Vector, Vector, Dune::OwnerOverlapCopyCommunication< int, int > > >, std::unique_ptr< typename EnergySolverSelector< Matrix, Vector >::type > > createParallelFlexibleSolver(const Grid &, const Matrix &, const PropertyTree &)
Definition: GenericTemperatureModel_impl.hpp:75
Definition: GenericTemperatureModel_impl.hpp:65
Dune::OverlappingSchwarzOperator< M, V, V, Comm > EnergyOperator
Definition: GenericTemperatureModel_impl.hpp:67
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: GenericTemperatureModel_impl.hpp:66