Go to the documentation of this file.
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>
39#include <opm/grid/CpGrid.hpp>
41#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
42#include <opm/input/eclipse/EclipseState/Tables/TracerVdTable.hpp>
43#include <opm/input/eclipse/Schedule/Well/Well.hpp>
44#include <opm/input/eclipse/Schedule/Well/WellTracerProperties.hpp>
53#include <fmt/format.h>
65template< class M, class V>
68 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
73template< class Vector, class Gr id, class Matrix>
74std::tuple<std::unique_ptr<Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
75 Dune::OwnerOverlapCopyCommunication<int,int>>>,
76 std::unique_ptr<typename TracerSolverSelector<Matrix,Vector>::type>>
79 OPM_THROW(std::logic_error, "Grid not supported for parallel Tracers.");
80 return { nullptr, nullptr};
83template< class Vector, class Matrix>
84std::tuple<std::unique_ptr<Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
85 Dune::OwnerOverlapCopyCommunication<int,int>>>,
86 std::unique_ptr<typename TracerSolverSelector<Matrix,Vector>::type>>
89 using TracerOperator = Dune::OverlappingSchwarzOperator<Matrix,Vector,Vector,
90 Dune::OwnerOverlapCopyCommunication<int,int>>;
92 const auto& cellComm = grid.cellCommunication();
93 auto op = std::make_unique<TracerOperator>(M, cellComm);
94 auto dummyWeights = [](){ return Vector();};
95 return {std::move(op), std::make_unique<TracerSolver>(*op, cellComm, prm, dummyWeights, 0)};
99template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
102 const EclipseState& eclState,
104 const DofMapper& dofMapper,
105 const std::function<std::array<double,dimWorld>( int)> centroids)
106 : gridView_(gridView)
107 , eclState_(eclState)
108 , cartMapper_(cartMapper)
109 , dofMapper_(dofMapper)
110 , centroids_(centroids)
114template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
118 if (freeTracerConcentration_.empty())
121 return freeTracerConcentration_[tracerIdx][globalDofIdx];
124template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
128 if (solTracerConcentration_.empty())
131 return solTracerConcentration_[tracerIdx][globalDofIdx];
134template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
138 this->freeTracerConcentration_[tracerIdx][globalDofIdx] = value;
139 this->tracerConcentration_[tracerIdx][globalDofIdx][0] = value;
142template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
146 this->solTracerConcentration_[tracerIdx][globalDofIdx] = value;
147 this->tracerConcentration_[tracerIdx][globalDofIdx][1] = value;
150template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
154 this->enableSolTracers_[tracerIdx] = enableSolTracer;
157template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
161 return this->eclState_.tracer().size();
164template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
166fname( int tracerIdx) const
168 return this->eclState_.tracer()[tracerIdx].fname();
171template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
173sname( int tracerIdx) const
175 return this->eclState_.tracer()[tracerIdx].sname();
178template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
182 return this->eclState_.tracer()[tracerIdx].wellfname();
185template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
189 return this->eclState_.tracer()[tracerIdx].wellsname();
192template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
194phase( int tracerIdx) const
196 return this->eclState_.tracer()[tracerIdx].phase;
199template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
203 return this->enableSolTracers_;
206template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
210 return eclWell.getTracerProperties().getConcentration(name);
213template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
215name( int tracerIdx) const
217 return this->eclState_.tracer()[tracerIdx].name;
220template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
222doInit( bool rst, std::size_t numGridDof,
223 std::size_t gasPhaseIdx, std::size_t oilPhaseIdx, std::size_t waterPhaseIdx)
225 const auto& tracers = eclState_.tracer();
227 if (tracers.size() == 0)
231 const std::size_t numTracers = tracers.size();
232 enableSolTracers_.resize(numTracers);
233 tracerConcentration_.resize(numTracers);
234 freeTracerConcentration_.resize(numTracers);
235 solTracerConcentration_.resize(numTracers);
238 tracerPhaseIdx_.resize(numTracers);
239 for (std::size_t tracerIdx = 0; tracerIdx < numTracers; tracerIdx++) {
240 const auto& tracer = tracers[tracerIdx];
242 if (tracer.phase == Phase::WATER)
243 tracerPhaseIdx_[tracerIdx] = waterPhaseIdx;
244 else if (tracer.phase == Phase::OIL)
245 tracerPhaseIdx_[tracerIdx] = oilPhaseIdx;
246 else if (tracer.phase == Phase::GAS)
247 tracerPhaseIdx_[tracerIdx] = gasPhaseIdx;
249 tracerConcentration_[tracerIdx].resize(numGridDof);
250 freeTracerConcentration_[tracerIdx].resize(numGridDof);
251 solTracerConcentration_[tracerIdx].resize(numGridDof);
258 if (tracer.free_concentration.has_value()){
259 const auto& free_concentration = tracer.free_concentration.value();
260 int tblkDatasize = free_concentration.size();
261 if (tblkDatasize < cartMapper_.cartesianSize()){
262 throw std::runtime_error( "Wrong size of TBLKF for" + tracer.name);
264 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
265 int cartDofIdx = cartMapper_.cartesianIndex(globalDofIdx);
266 tracerConcentration_[tracerIdx][globalDofIdx][0] = free_concentration[cartDofIdx];
267 freeTracerConcentration_[tracerIdx][globalDofIdx] = free_concentration[cartDofIdx];
271 else if (tracer.free_tvdp.has_value()) {
272 const auto& free_tvdp = tracer.free_tvdp.value();
273 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
274 tracerConcentration_[tracerIdx][globalDofIdx][0] =
275 free_tvdp.evaluate( "TRACER_CONCENTRATION",
276 centroids_(globalDofIdx)[2]);
277 freeTracerConcentration_[tracerIdx][globalDofIdx] =
278 free_tvdp.evaluate( "TRACER_CONCENTRATION",
279 centroids_(globalDofIdx)[2]);
283 OpmLog::warning(fmt::format( "No TBLKF or TVDPF given for free tracer {}. "
284 "Initial values set to zero. ", tracer.name));
285 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
286 tracerConcentration_[tracerIdx][globalDofIdx][0] = 0.0;
287 freeTracerConcentration_[tracerIdx][globalDofIdx] = 0.0;
292 if (tracer.phase != Phase::WATER &&
293 ((tracer.phase == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
294 (tracer.phase == Phase::OIL && FluidSystem::enableVaporizedOil()))) {
296 if (tracer.solution_concentration.has_value()){
297 enableSolTracers_[tracerIdx] = true;
298 const auto& solution_concentration = tracer.solution_concentration.value();
299 int tblkDatasize = solution_concentration.size();
300 if (tblkDatasize < cartMapper_.cartesianSize()){
301 throw std::runtime_error( "Wrong size of TBLKS for" + tracer.name);
303 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
304 int cartDofIdx = cartMapper_.cartesianIndex(globalDofIdx);
305 tracerConcentration_[tracerIdx][globalDofIdx][1] = solution_concentration[cartDofIdx];
306 solTracerConcentration_[tracerIdx][globalDofIdx] = solution_concentration[cartDofIdx];
310 else if (tracer.solution_tvdp.has_value()) {
311 enableSolTracers_[tracerIdx] = true;
312 const auto& solution_tvdp = tracer.solution_tvdp.value();
313 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
314 tracerConcentration_[tracerIdx][globalDofIdx][1] =
315 solution_tvdp.evaluate( "TRACER_CONCENTRATION",
316 centroids_(globalDofIdx)[2]);
317 solTracerConcentration_[tracerIdx][globalDofIdx] =
318 solution_tvdp.evaluate( "TRACER_CONCENTRATION",
319 centroids_(globalDofIdx)[2]);
324 enableSolTracers_[tracerIdx] = false;
325 OpmLog::warning(fmt::format( "No TBLKS or TVDPS given for solution tracer {}. "
326 "Initial values set to zero. ", tracer.name));
327 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
328 tracerConcentration_[tracerIdx][globalDofIdx][1] = 0.0;
329 solTracerConcentration_[tracerIdx][globalDofIdx] = 0.0;
335 enableSolTracers_[tracerIdx] = false;
336 for (std::size_t globalDofIdx = 0; globalDofIdx < numGridDof; ++globalDofIdx) {
337 tracerConcentration_[tracerIdx][globalDofIdx][1] = 0.0;
338 solTracerConcentration_[tracerIdx][globalDofIdx] = 0.0;
344 tracerMatrix_ = std::make_unique<TracerMatrix>(numGridDof, numGridDof, TracerMatrix::random);
347 using NeighborSet = std::set<unsigned>;
348 std::vector<NeighborSet> neighbors(numGridDof);
350 Stencil stencil(gridView_, dofMapper_);
351 for ( const auto& elem : elements(gridView_)) {
352 stencil.update(elem);
354 for ( unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
355 unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
357 for ( unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
358 unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
359 neighbors[myIdx].insert(neighborIdx);
365 for ( unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
366 tracerMatrix_->setrowsize(dofIdx, neighbors[dofIdx].size());
368 tracerMatrix_->endrowsizes();
373 for ( unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
374 typename NeighborSet::iterator nIt = neighbors[dofIdx].begin();
375 typename NeighborSet::iterator nEndIt = neighbors[dofIdx].end();
376 for (; nIt != nEndIt; ++nIt) {
377 tracerMatrix_->addindex(dofIdx, *nIt);
380 tracerMatrix_->endindices();
383template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
388 Scalar tolerance = 1e-2;
393 prm. put( "maxiter", maxIter);
394 prm. put( "tol", tolerance);
395 prm. put( "verbosity", verbosity);
396 prm. put( "solver", std::string( "bicgstab"));
397 prm. put( "preconditioner.type", std::string( "ParOverILU0"));
400 if(gridView_.grid().comm().size() > 1)
402 auto [tracerOperator, solver] =
403 createParallelFlexibleSolver<TracerVector>(gridView_.grid(), M, prm);
404 (void) tracerOperator;
407 solver->apply(x, b, result);
410 return result.converged;
415 using TracerSolver = Dune::BiCGSTABSolver<TracerVector>;
416 using TracerOperator = Dune::MatrixAdapter<TracerMatrix,TracerVector,TracerVector>;
417 using TracerScalarProduct = Dune::SeqScalarProduct<TracerVector>;
418 using TracerPreconditioner = Dune::SeqILU< TracerMatrix,TracerVector,TracerVector>;
420 TracerOperator tracerOperator(M);
421 TracerScalarProduct tracerScalarProduct;
422 TracerPreconditioner tracerPreconditioner(M, 0, 1);
424 TracerSolver solver (tracerOperator, tracerScalarProduct,
425 tracerPreconditioner, tolerance, maxIter,
429 solver.apply(x, b, result);
432 return result.converged;
438template< class Gr id, class Gr idView, class DofMapper, class Stencil, class Flu idSystem, class Scalar>
442 Scalar tolerance = 1e-2;
447 prm. put( "maxiter", maxIter);
448 prm. put( "tol", tolerance);
449 prm. put( "verbosity", verbosity);
450 prm. put( "solver", std::string( "bicgstab"));
451 prm. put( "preconditioner.type", std::string( "ParOverILU0"));
454 if(gridView_.grid().comm().size() > 1)
456 auto [tracerOperator, solver] =
457 createParallelFlexibleSolver<TracerVector>(gridView_.grid(), M, prm);
458 (void) tracerOperator;
459 bool converged = true;
460 for (std::size_t nrhs = 0; nrhs < b.size(); ++nrhs) {
463 solver->apply(x[nrhs], b[nrhs], result);
464 converged = (converged && result.converged);
471 using TracerSolver = Dune::BiCGSTABSolver<TracerVector>;
472 using TracerOperator = Dune::MatrixAdapter<TracerMatrix,TracerVector,TracerVector>;
473 using TracerScalarProduct = Dune::SeqScalarProduct<TracerVector>;
474 using TracerPreconditioner = Dune::SeqILU< TracerMatrix,TracerVector,TracerVector>;
476 TracerOperator tracerOperator(M);
477 TracerScalarProduct tracerScalarProduct;
478 TracerPreconditioner tracerPreconditioner(M, 0, 1);
480 TracerSolver solver (tracerOperator, tracerScalarProduct,
481 tracerPreconditioner, tolerance, maxIter,
484 bool converged = true;
485 for (std::size_t nrhs = 0; nrhs < b.size(); ++nrhs) {
488 solver.apply(x[nrhs], b[nrhs], result);
489 converged = (converged && result.converged);
Definition: CollectDataOnIORank.hpp:49
Definition: FlexibleSolver.hpp:45
void setEnableSolTracers(int tracerIdx, bool enableSolTracer) Definition: GenericTracerModel_impl.hpp:152
int numTracers() const Return the number of tracers considered by the tracerModel. Definition: GenericTracerModel_impl.hpp:159
Scalar currentConcentration_(const Well &eclWell, const std::string &name) const Definition: GenericTracerModel_impl.hpp:208
void doInit(bool rst, std::size_t numGridDof, std::size_t gasPhaseIdx, std::size_t oilPhaseIdx, std::size_t waterPhaseIdx) Initialize all internal data structures needed by the tracer module. Definition: GenericTracerModel_impl.hpp:222
std::string sname(int tracerIdx) const Definition: GenericTracerModel_impl.hpp:173
void setSolTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value) Definition: GenericTracerModel_impl.hpp:144
std::string wellfname(int tracerIdx) const Definition: GenericTracerModel_impl.hpp:180
GenericTracerModel(const GridView &gridView, const EclipseState &eclState, const CartesianIndexMapper &cartMapper, const DofMapper &dofMapper, const std::function< std::array< double, dimWorld >(int)> centroids) Definition: GenericTracerModel_impl.hpp:101
const std::string & name(int tracerIdx) const Return the tracer name. Definition: GenericTracerModel_impl.hpp:215
const std::vector< bool > & enableSolTracers() const Definition: GenericTracerModel_impl.hpp:201
Phase phase(int tracerIdx) const Definition: GenericTracerModel_impl.hpp:194
std::string wellsname(int tracerIdx) const Definition: GenericTracerModel_impl.hpp:187
std::string fname(int tracerIdx) const Definition: GenericTracerModel_impl.hpp:166
bool linearSolve_(const TracerMatrix &M, TracerVector &x, TracerVector &b) Definition: GenericTracerModel_impl.hpp:385
Dune::BlockVector< Dune::FieldVector< Scalar, 2 > > TracerVector Definition: GenericTracerModel.hpp:59
void setFreeTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value) Definition: GenericTracerModel_impl.hpp:136
Scalar solTracerConcentration(int tracerIdx, int globalDofIdx) const Definition: GenericTracerModel_impl.hpp:126
bool linearSolveBatchwise_(const TracerMatrix &M, std::vector< TracerVector > &x, std::vector< TracerVector > &b) Definition: GenericTracerModel_impl.hpp:440
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 2, 2 > > TracerMatrix Definition: GenericTracerModel.hpp:58
Scalar freeTracerConcentration(int tracerIdx, int globalDofIdx) const Return the tracer concentration for tracer index and global DofIdx. Definition: GenericTracerModel_impl.hpp:116
Definition: PropertyTree.hpp:37
void put(const std::string &key, const T &data)
Definition: blackoilboundaryratevector.hh:37
Dune::InverseOperatorResult InverseOperatorResult Definition: BdaBridge.hpp:32
std::tuple< std::unique_ptr< Dune::OverlappingSchwarzOperator< Matrix, Vector, Vector, Dune::OwnerOverlapCopyCommunication< int, int > > >, std::unique_ptr< typename TracerSolverSelector< Matrix, Vector >::type > > createParallelFlexibleSolver(const Grid &, const Matrix &, const PropertyTree &) Definition: GenericTracerModel_impl.hpp:77
Definition: GenericTracerModel_impl.hpp:67
Dune::FlexibleSolver< TracerOperator > type Definition: GenericTracerModel_impl.hpp:70
Dune::OverlappingSchwarzOperator< M, V, V, Comm > TracerOperator Definition: GenericTracerModel_impl.hpp:69
Dune::OwnerOverlapCopyCommunication< int, int > Comm Definition: GenericTracerModel_impl.hpp:68
|