28#ifndef OPM_GENERIC_TRACER_MODEL_IMPL_HPP
29#define OPM_GENERIC_TRACER_MODEL_IMPL_HPP
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>
37#include <opm/common/OpmLog/OpmLog.hpp>
38#include <opm/common/TimingMacros.hpp>
40#include <opm/grid/CpGrid.hpp>
42#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
43#include <opm/input/eclipse/EclipseState/Tables/TracerVdTable.hpp>
44#include <opm/input/eclipse/Schedule/Well/Well.hpp>
45#include <opm/input/eclipse/Schedule/Well/WellTracerProperties.hpp>
50#include <opm/simulators/linalg/ilufirstelement.hh>
51#include <opm/simulators/linalg/PropertyTree.hpp>
52#include <opm/simulators/linalg/FlexibleSolver.hpp>
54#include <fmt/format.h>
66template<
class M,
class V>
67struct TracerSolverSelector
69 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
70 using TracerOperator = Dune::OverlappingSchwarzOperator<M, V, V, Comm>;
71 using type = Dune::FlexibleSolver<TracerOperator>;
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 TracerSolverSelector<Matrix,Vector>::type>>
78createParallelFlexibleSolver(
const Grid&,
const Matrix&,
const PropertyTree&)
80 OPM_THROW(std::logic_error,
"Grid not supported for parallel Tracers.");
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 TracerSolverSelector<Matrix,Vector>::type>>
88createParallelFlexibleSolver(
const Dune::CpGrid& grid,
const Matrix& M,
const PropertyTree& prm)
90 using TracerOperator = Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
91 Dune::OwnerOverlapCopyCommunication<int,int>>;
92 using TracerSolver = Dune::FlexibleSolver<TracerOperator>;
93 const auto& cellComm = grid.cellCommunication();
94 auto op = std::make_unique<TracerOperator>(M, cellComm);
95 auto dummyWeights = [](){
return Vector();};
96 return {std::move(op), std::make_unique<TracerSolver>(*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 const std::function<std::array<double,dimWorld>(
int)> centroids)
107 : gridView_(gridView)
108 , eclState_(eclState)
109 , cartMapper_(cartMapper)
110 , dofMapper_(dofMapper)
111 , centroids_(centroids)
115template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
116Scalar GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
117freeTracerConcentration(
int tracerIdx,
int globalDofIdx)
const
119 if (freeTracerConcentration_.empty()) {
123 return freeTracerConcentration_[tracerIdx][globalDofIdx];
126template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
130 if (solTracerConcentration_.empty()) {
134 return solTracerConcentration_[tracerIdx][globalDofIdx];
137template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
138void GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
139setFreeTracerConcentration(
int tracerIdx,
int globalDofIdx, Scalar value)
141 this->freeTracerConcentration_[tracerIdx][globalDofIdx] = value;
142 this->tracerConcentration_[tracerIdx][globalDofIdx][0] = value;
145template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
146void GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
147setSolTracerConcentration(
int tracerIdx,
int globalDofIdx, Scalar value)
149 this->solTracerConcentration_[tracerIdx][globalDofIdx] = value;
150 this->tracerConcentration_[tracerIdx][globalDofIdx][1] = value;
153template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
154void GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
155setEnableSolTracers(
int tracerIdx,
bool enableSolTracer)
157 this->enableSolTracers_[tracerIdx] = enableSolTracer;
160template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
161int GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
164 return this->eclState_.tracer().size();
167template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
169fname(
int tracerIdx)
const
171 return this->eclState_.tracer()[tracerIdx].fname();
174template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
175std::string GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
176sname(
int tracerIdx)
const
178 return this->eclState_.tracer()[tracerIdx].sname();
181template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
182std::string GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
183wellfname(
int tracerIdx)
const
185 return this->eclState_.tracer()[tracerIdx].wellfname();
188template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
189std::string GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
190wellsname(
int tracerIdx)
const
192 return this->eclState_.tracer()[tracerIdx].wellsname();
195template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
196Phase GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
197phase(
int tracerIdx)
const
199 return this->eclState_.tracer()[tracerIdx].phase;
202template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
203const std::vector<bool>& GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
204enableSolTracers()
const
206 return this->enableSolTracers_;
209template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
210Scalar GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
211currentConcentration_(
const Well& eclWell,
const std::string& trName,
const SummaryState& summaryState)
const
213 return eclWell.getTracerProperties().getConcentration(WellTracerProperties::Well { eclWell.name() },
214 WellTracerProperties::Tracer { trName },
218template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
219const std::string& GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
220name(
int tracerIdx)
const
222 return this->eclState_.tracer()[tracerIdx].name;
225template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
226void GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
227doInit(
bool rst, std::size_t numGridDof,
228 std::size_t gasPhaseIdx, std::size_t oilPhaseIdx, std::size_t waterPhaseIdx)
230 const auto& tracers = eclState_.tracer();
232 if (tracers.size() == 0) {
237 const std::size_t
numTracers = tracers.size();
245 for (std::size_t tracerIdx = 0; tracerIdx <
numTracers; tracerIdx++) {
246 const auto& tracer = tracers[tracerIdx];
248 if (tracer.phase == Phase::WATER)
249 tracerPhaseIdx_[tracerIdx] = waterPhaseIdx;
250 else if (tracer.phase == Phase::OIL)
251 tracerPhaseIdx_[tracerIdx] = oilPhaseIdx;
252 else if (tracer.phase == Phase::GAS)
253 tracerPhaseIdx_[tracerIdx] = gasPhaseIdx;
255 tracerConcentration_[tracerIdx].resize(numGridDof);
256 freeTracerConcentration_[tracerIdx].resize(numGridDof);
257 solTracerConcentration_[tracerIdx].resize(numGridDof);
264 if (tracer.free_concentration.has_value()){
265 const auto& free_concentration = tracer.free_concentration.value();
266 int tblkDatasize = free_concentration.size();
267 if (tblkDatasize < cartMapper_.cartesianSize()){
268 throw std::runtime_error(
"Wrong size of TBLKF for" + tracer.name);
270 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
271 int cartDofIdx = cartMapper_.cartesianIndex(globalDofIdx);
272 tracerConcentration_[tracerIdx][globalDofIdx][0] = free_concentration[cartDofIdx];
273 freeTracerConcentration_[tracerIdx][globalDofIdx] = free_concentration[cartDofIdx];
277 else if (tracer.free_tvdp.has_value()) {
278 const auto& free_tvdp = tracer.free_tvdp.value();
279 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
280 tracerConcentration_[tracerIdx][globalDofIdx][0] =
281 free_tvdp.evaluate(
"TRACER_CONCENTRATION",
283 freeTracerConcentration_[tracerIdx][globalDofIdx] =
284 free_tvdp.evaluate(
"TRACER_CONCENTRATION",
289 OpmLog::warning(fmt::format(
"No TBLKF or TVDPF given for free tracer {}. "
290 "Initial values set to zero. ", tracer.name));
291 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
292 tracerConcentration_[tracerIdx][globalDofIdx][0] = 0.0;
293 freeTracerConcentration_[tracerIdx][globalDofIdx] = 0.0;
298 if (tracer.phase != Phase::WATER &&
299 ((tracer.phase == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
300 (tracer.phase == Phase::OIL && FluidSystem::enableVaporizedOil()))) {
302 if (tracer.solution_concentration.has_value()){
303 enableSolTracers_[tracerIdx] =
true;
304 const auto& solution_concentration = tracer.solution_concentration.value();
305 int tblkDatasize = solution_concentration.size();
306 if (tblkDatasize < cartMapper_.cartesianSize()){
307 throw std::runtime_error(
"Wrong size of TBLKS for" + tracer.name);
309 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
310 int cartDofIdx = cartMapper_.cartesianIndex(globalDofIdx);
311 tracerConcentration_[tracerIdx][globalDofIdx][1] = solution_concentration[cartDofIdx];
312 solTracerConcentration_[tracerIdx][globalDofIdx] = solution_concentration[cartDofIdx];
316 else if (tracer.solution_tvdp.has_value()) {
317 enableSolTracers_[tracerIdx] =
true;
318 const auto& solution_tvdp = tracer.solution_tvdp.value();
319 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
320 tracerConcentration_[tracerIdx][globalDofIdx][1] =
321 solution_tvdp.evaluate(
"TRACER_CONCENTRATION",
323 solTracerConcentration_[tracerIdx][globalDofIdx] =
324 solution_tvdp.evaluate(
"TRACER_CONCENTRATION",
330 enableSolTracers_[tracerIdx] =
false;
331 OpmLog::warning(fmt::format(
"No TBLKS or TVDPS given for solution tracer {}. "
332 "Initial values set to zero. ", tracer.name));
333 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
334 tracerConcentration_[tracerIdx][globalDofIdx][1] = 0.0;
335 solTracerConcentration_[tracerIdx][globalDofIdx] = 0.0;
341 enableSolTracers_[tracerIdx] =
false;
342 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
343 tracerConcentration_[tracerIdx][globalDofIdx][1] = 0.0;
344 solTracerConcentration_[tracerIdx][globalDofIdx] = 0.0;
350 tracerMatrix_ = std::make_unique<TracerMatrix>(numGridDof, numGridDof, TracerMatrix::random);
353 using NeighborSet = std::set<unsigned>;
354 std::vector<NeighborSet> neighbors(numGridDof);
356 Stencil stencil(gridView_, dofMapper_);
357 for (
const auto& elem : elements(gridView_)) {
358 stencil.update(elem);
360 for (
unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
361 unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
363 for (
unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
364 unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
365 neighbors[myIdx].insert(neighborIdx);
371 for (
unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
372 tracerMatrix_->setrowsize(dofIdx, neighbors[dofIdx].size());
374 tracerMatrix_->endrowsizes();
379 for (
unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
380 for (
const auto& index : neighbors[dofIdx]) {
381 tracerMatrix_->addindex(dofIdx, index);
384 tracerMatrix_->endindices();
387template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
389linearSolve_(
const TracerMatrix& M, TracerVector& x, TracerVector& b)
392 Scalar tolerance = 1e-2;
397 prm.
put(
"maxiter", maxIter);
398 prm.
put(
"tol", tolerance);
399 prm.
put(
"verbosity", verbosity);
400 prm.
put(
"solver", std::string(
"bicgstab"));
401 prm.
put(
"preconditioner.type", std::string(
"paroverilu0"));
404 if(gridView_.grid().comm().size() > 1)
406 auto [tracerOperator, solver] =
407 createParallelFlexibleSolver<TracerVector>(gridView_.grid(), M, prm);
408 (void) tracerOperator;
410 Dune::InverseOperatorResult result;
411 solver->apply(x, b, result);
414 return result.converged;
419 using TracerSolver = Dune::BiCGSTABSolver<TracerVector>;
420 using TracerOperator = Dune::MatrixAdapter<TracerMatrix,TracerVector,TracerVector>;
421 using TracerScalarProduct = Dune::SeqScalarProduct<TracerVector>;
422 using TracerPreconditioner = Dune::SeqILU< TracerMatrix,TracerVector,TracerVector>;
424 TracerOperator tracerOperator(M);
425 TracerScalarProduct tracerScalarProduct;
426 TracerPreconditioner tracerPreconditioner(M, 0, 1);
428 TracerSolver solver (tracerOperator, tracerScalarProduct,
429 tracerPreconditioner, tolerance, maxIter,
432 Dune::InverseOperatorResult result;
433 solver.apply(x, b, result);
436 return result.converged;
442template<
class Gr
id,
class Gr
idView,
class DofMapper,
class Stencil,
class Flu
idSystem,
class Scalar>
443bool GenericTracerModel<Grid,GridView,DofMapper,Stencil,FluidSystem,Scalar>::
444linearSolveBatchwise_(
const TracerMatrix& M, std::vector<TracerVector>& x, std::vector<TracerVector>& b)
446 OPM_TIMEBLOCK(tracerSolve);
447 const Scalar tolerance = 1e-2;
448 const int maxIter = 100;
449 const int verbosity = 0;
452 if (gridView_.grid().comm().size() > 1)
455 prm.
put(
"maxiter", maxIter);
456 prm.put(
"tol", tolerance);
457 prm.put(
"verbosity", verbosity);
458 prm.put(
"solver", std::string(
"bicgstab"));
459 prm.put(
"preconditioner.type", std::string(
"paroverilu0"));
460 auto [tracerOperator, solver] =
461 createParallelFlexibleSolver<TracerVector>(gridView_.grid(), M, prm);
462 (void) tracerOperator;
463 bool converged =
true;
464 for (std::size_t nrhs = 0; nrhs < b.size(); ++nrhs) {
466 Dune::InverseOperatorResult result;
467 solver->apply(x[nrhs], b[nrhs], result);
468 converged = (converged && result.converged);
475 using TracerSolver = Dune::BiCGSTABSolver<TracerVector>;
476 using TracerOperator = Dune::MatrixAdapter<TracerMatrix,TracerVector,TracerVector>;
477 using TracerScalarProduct = Dune::SeqScalarProduct<TracerVector>;
478 using TracerPreconditioner = Dune::SeqILU<TracerMatrix,TracerVector,TracerVector>;
480 if (std::ranges::all_of(b, [](
const auto& v) {
return v.infinity_norm() == 0.0; }))
485 TracerOperator tracerOperator(M);
486 TracerScalarProduct tracerScalarProduct;
487 TracerPreconditioner tracerPreconditioner(M, 0, 1);
489 TracerSolver solver (tracerOperator, tracerScalarProduct,
490 tracerPreconditioner, tolerance, maxIter,
493 bool converged =
true;
494 for (std::size_t nrhs = 0; nrhs < b.size(); ++nrhs) {
496 Dune::InverseOperatorResult result;
497 solver.apply(x[nrhs], b[nrhs], result);
498 converged = (converged && result.converged);
A class which handles tracers as specified in by ECL.
Definition GenericTracerModel.hpp:56
int numTracers() const
Return the number of tracers considered by the tracerModel.
Definition GenericTracerModel_impl.hpp:162
std::function< std::array< double, dimWorld >(int)> centroids_
Function returning the cell centers.
Definition GenericTracerModel.hpp:171
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