28#ifndef OPM_GENERIC_TEMPERATURE_MODEL_IMPL_HPP
29#define OPM_GENERIC_TEMPERATURE_MODEL_IMPL_HPP
32#ifndef OPM_GENERIC_TEMPERATURE_MODEL_HPP
37#include <dune/istl/operators.hh>
38#include <dune/istl/solvers.hh>
39#include <dune/istl/schwarz.hh>
40#include <dune/istl/preconditioners.hh>
41#include <dune/istl/schwarz.hh>
43#include <opm/common/OpmLog/OpmLog.hpp>
45#include <opm/grid/CpGrid.hpp>
47#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
48#include <opm/input/eclipse/Schedule/Well/Well.hpp>
52#include <opm/simulators/linalg/ilufirstelement.hh>
53#include <opm/simulators/linalg/PropertyTree.hpp>
54#include <opm/simulators/linalg/FlexibleSolver.hpp>
61#include <fmt/format.h>
66template<
class M,
class V>
67struct EnergySolverSelector
69 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
70 using EnergyOperator = Dune::OverlappingSchwarzOperator<M, V, V, Comm>;
71 using type = Dune::FlexibleSolver<EnergyOperator>;
74template<
class Vector,
class Gr
id,
class Matrix>
75std::tuple<std::unique_ptr<Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
76 Dune::OwnerOverlapCopyCommunication<int,int>>>,
77 std::unique_ptr<typename EnergySolverSelector<Matrix,Vector>::type>>
78createParallelFlexibleSolver(
const Grid&,
const Matrix&,
const PropertyTree&)
80 OPM_THROW(std::logic_error,
"Grid not supported for parallel Temperatures.");
81 return {
nullptr,
nullptr};
84template<
class Vector,
class Matrix>
85std::tuple<std::unique_ptr<Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
86 Dune::OwnerOverlapCopyCommunication<int,int>>>,
87 std::unique_ptr<typename EnergySolverSelector<Matrix,Vector>::type>>
88createParallelFlexibleSolver(
const Dune::CpGrid& grid,
const Matrix& M,
const PropertyTree& prm)
90 using EnergyOperator = Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
91 Dune::OwnerOverlapCopyCommunication<int,int>>;
92 using EnergySolver = Dune::FlexibleSolver<EnergyOperator>;
93 const auto& cellComm = grid.cellCommunication();
94 auto op = std::make_unique<EnergyOperator>(M, cellComm);
95 auto dummyWeights = [](){
return Vector();};
96 return {std::move(op), std::make_unique<EnergySolver>(*op, cellComm, prm, dummyWeights, 0)};
100template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
103 const EclipseState& eclState,
104 const CartesianIndexMapper& cartMapper,
105 const DofMapper& dofMapper)
106 : gridView_(gridView)
107 , eclState_(eclState)
108 , cartMapper_(cartMapper)
109 , dofMapper_(dofMapper)
113template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
114void GenericTemperatureModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
115doInit(std::size_t numGridDof)
117 doTemp_ = eclState_.getSimulationConfig().isTemp();
118 temperature_.resize(numGridDof);
123 energyVector_.resize(numGridDof);
127template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
131 Scalar tolerance = 1e-2;
135 prm.
put(
"maxiter", maxIter);
136 prm.
put(
"tol", tolerance);
137 prm.
put(
"verbosity", verbosity);
138 prm.
put(
"solver", std::string(
"bicgstab"));
139 prm.
put(
"preconditioner.type", std::string(
"ParOverILU0"));
142 if (gridView_.grid().comm().size() > 1)
144 auto [energyOperator, solver] =
145 createParallelFlexibleSolver<EnergyVector>(gridView_.grid(), matrix, prm);
146 op_ = std::move(energyOperator);
147 pre_ = &solver->preconditioner();
148 linear_solver_ = std::move(solver);
153 using SeqOperator = Dune::MatrixAdapter<EnergyMatrix, EnergyVector, EnergyVector>;
155 auto seqOp = std::make_unique<SeqOperator>(matrix);
156 auto dummyWeights = [](){
return EnergyVector(); };
157 auto seqSolver = std::make_unique<SeqSolver>(*seqOp, prm, dummyWeights, 0);
158 op_ = std::move(seqOp);
159 pre_ = &seqSolver->preconditioner();
160 linear_solver_ = std::move(seqSolver);
164template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
165bool GenericTemperatureModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
166linearSolve_(
const EnergyMatrix& M, EnergyVector& x, EnergyVector& b)
168 if (!linear_solver_) {
169 setupLinearSolver(M);
174 Dune::InverseOperatorResult result;
175 linear_solver_->apply(x, b, result);
176 return result.converged;
A class which handles sequential implicit solution of the energy equation as specified in by TEMP.
A solver class that encapsulates all needed objects for a linear solver (operator,...
Definition FlexibleSolver.hpp:45
Definition GenericTemperatureModel.hpp:53
Hierarchical collection of key/value pairs.
Definition PropertyTree.hpp:39
void put(const std::string &key, const T &data)
Insert key/value pair into property tree.
Definition PropertyTree.cpp:71
Represents the stencil (finite volume geometry) of a single element in the ECFV discretization.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187