GenericTracerModel_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_TRACER_MODEL_IMPL_HPP
29#define OPM_GENERIC_TRACER_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#include <opm/common/TimingMacros.hpp>
39
40#include <opm/grid/CpGrid.hpp>
41
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>
46
48
53
54#include <fmt/format.h>
55
56#include <array>
57#include <functional>
58#include <memory>
59#include <set>
60#include <stdexcept>
61#include <string>
62
63namespace Opm {
64
65#if HAVE_MPI
66template<class M, class V>
68{
69 using Comm = Dune::OwnerOverlapCopyCommunication<int, int>;
70 using TracerOperator = Dune::OverlappingSchwarzOperator<M, V, V, Comm>;
72};
74template<class Vector, class Grid, 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&)
79{
80 OPM_THROW(std::logic_error, "Grid not supported for parallel Tracers.");
81 return {nullptr, nullptr};
82}
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)
89{
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)};
97}
98#endif
99
100template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
102GenericTracerModel(const GridView& gridView,
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)
112{
113}
114
115template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
117freeTracerConcentration(int tracerIdx, int globalDofIdx) const
118{
119 if (freeTracerConcentration_.empty()) {
120 return 0.0;
121 }
122
123 return freeTracerConcentration_[tracerIdx][globalDofIdx];
124}
125
126template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
128solTracerConcentration(int tracerIdx, int globalDofIdx) const
130 if (solTracerConcentration_.empty()) {
131 return 0.0;
132 }
133
134 return solTracerConcentration_[tracerIdx][globalDofIdx];
136
137template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
139setFreeTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
140{
141 this->freeTracerConcentration_[tracerIdx][globalDofIdx] = value;
142 this->tracerConcentration_[tracerIdx][globalDofIdx][0] = value;
143}
144
145template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
147setSolTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
148{
149 this->solTracerConcentration_[tracerIdx][globalDofIdx] = value;
150 this->tracerConcentration_[tracerIdx][globalDofIdx][1] = value;
151}
152
153template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
155setEnableSolTracers(int tracerIdx, bool enableSolTracer)
156{
157 this->enableSolTracers_[tracerIdx] = enableSolTracer;
158}
159
160template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
162numTracers() const
163{
164 return this->eclState_.tracer().size();
165}
166
167template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
169fname(int tracerIdx) const
170{
171 return this->eclState_.tracer()[tracerIdx].fname();
172}
173
174template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
176sname(int tracerIdx) const
177{
178 return this->eclState_.tracer()[tracerIdx].sname();
179}
180
181template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
183wellfname(int tracerIdx) const
184{
185 return this->eclState_.tracer()[tracerIdx].wellfname();
186}
187
188template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
190wellsname(int tracerIdx) const
191{
192 return this->eclState_.tracer()[tracerIdx].wellsname();
193}
194
195template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
197phase(int tracerIdx) const
198{
199 return this->eclState_.tracer()[tracerIdx].phase;
200}
201
202template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
204enableSolTracers() const
205{
206 return this->enableSolTracers_;
207}
208
209template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
211currentConcentration_(const Well& eclWell, const std::string& trName, const SummaryState& summaryState) const
212{
213 return eclWell.getTracerProperties().getConcentration(WellTracerProperties::Well { eclWell.name() },
214 WellTracerProperties::Tracer { trName },
215 summaryState);
216}
217
218template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
220name(int tracerIdx) const
221{
222 return this->eclState_.tracer()[tracerIdx].name;
223}
224
225template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
227doInit(bool rst, std::size_t numGridDof,
228 std::size_t gasPhaseIdx, std::size_t oilPhaseIdx, std::size_t waterPhaseIdx)
229{
230 const auto& tracers = eclState_.tracer();
231
232 if (tracers.size() == 0) {
233 return; // tracer treatment is supposed to be disabled
234 }
235
236 // retrieve the number of tracers from the deck
237 const std::size_t numTracers = tracers.size();
238 enableSolTracers_.resize(numTracers);
239 tracerConcentration_.resize(numTracers);
240 freeTracerConcentration_.resize(numTracers);
241 solTracerConcentration_.resize(numTracers);
242
243 // the phase where the tracer is
244 tracerPhaseIdx_.resize(numTracers);
245 for (std::size_t tracerIdx = 0; tracerIdx < numTracers; tracerIdx++) {
246 const auto& tracer = tracers[tracerIdx];
247
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;
254
255 tracerConcentration_[tracerIdx].resize(numGridDof);
256 freeTracerConcentration_[tracerIdx].resize(numGridDof);
257 solTracerConcentration_[tracerIdx].resize(numGridDof);
258
259 if (rst)
260 continue;
261
262
263 // TBLKF keyword
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);
269 }
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];
274 }
275 }
276 // TVDPF keyword
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",
282 centroids_(globalDofIdx)[2]);
283 freeTracerConcentration_[tracerIdx][globalDofIdx] =
284 free_tvdp.evaluate("TRACER_CONCENTRATION",
285 centroids_(globalDofIdx)[2]);
286 }
287 }
288 else {
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;
294 }
295 }
296
297 // Solution tracer initialization only needed for gas/oil tracers with DISGAS/VAPOIL active
298 if (tracer.phase != Phase::WATER &&
299 ((tracer.phase == Phase::GAS && FluidSystem::enableDissolvedGas()) ||
300 (tracer.phase == Phase::OIL && FluidSystem::enableVaporizedOil()))) {
301 // TBLKS keyword
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);
308 }
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];
313 }
314 }
315 // TVDPS keyword
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",
322 centroids_(globalDofIdx)[2]);
323 solTracerConcentration_[tracerIdx][globalDofIdx] =
324 solution_tvdp.evaluate("TRACER_CONCENTRATION",
325 centroids_(globalDofIdx)[2]);
326 }
327 }
328 else {
329 // No solution tracers, default to zero
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;
336 }
337 }
338 }
339 else {
340 // No solution tracers, default to zero
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;
345 }
346 }
347 }
348
349 // allocate matrix for storing the Jacobian of the tracer residual
350 tracerMatrix_ = std::make_unique<TracerMatrix>(numGridDof, numGridDof, TracerMatrix::random);
351
352 // find the sparsity pattern of the tracer matrix
353 using NeighborSet = std::set<unsigned>;
354 std::vector<NeighborSet> neighbors(numGridDof);
355
356 Stencil stencil(gridView_, dofMapper_);
357 for (const auto& elem : elements(gridView_)) {
358 stencil.update(elem);
359
360 for (unsigned primaryDofIdx = 0; primaryDofIdx < stencil.numPrimaryDof(); ++primaryDofIdx) {
361 unsigned myIdx = stencil.globalSpaceIndex(primaryDofIdx);
362
363 for (unsigned dofIdx = 0; dofIdx < stencil.numDof(); ++dofIdx) {
364 unsigned neighborIdx = stencil.globalSpaceIndex(dofIdx);
365 neighbors[myIdx].insert(neighborIdx);
366 }
367 }
368 }
369
370 // allocate space for the rows of the matrix
371 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
372 tracerMatrix_->setrowsize(dofIdx, neighbors[dofIdx].size());
373 }
374 tracerMatrix_->endrowsizes();
375
376 // fill the rows with indices. each degree of freedom talks to
377 // all of its neighbors. (it also talks to itself since
378 // degrees of freedom are sometimes quite egocentric.)
379 for (unsigned dofIdx = 0; dofIdx < numGridDof; ++ dofIdx) {
380 for (const auto& index : neighbors[dofIdx]) {
381 tracerMatrix_->addindex(dofIdx, index);
382 }
383 }
384 tracerMatrix_->endindices();
385}
386
387template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
390{
391 x = 0.0;
392 Scalar tolerance = 1e-2;
393 int maxIter = 100;
394
395 int verbosity = 0;
396 PropertyTree prm;
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"));
402
403#if HAVE_MPI
404 if(gridView_.grid().comm().size() > 1)
405 {
406 auto [tracerOperator, solver] =
407 createParallelFlexibleSolver<TracerVector>(gridView_.grid(), M, prm);
408 (void) tracerOperator;
409
411 solver->apply(x, b, result);
412
413 // return the result of the solver
414 return result.converged;
415 }
416 else
417 {
418#endif
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>;
423
424 TracerOperator tracerOperator(M);
425 TracerScalarProduct tracerScalarProduct;
426 TracerPreconditioner tracerPreconditioner(M, 0, 1); // results in ILU0
427
428 TracerSolver solver (tracerOperator, tracerScalarProduct,
429 tracerPreconditioner, tolerance, maxIter,
430 verbosity);
431
433 solver.apply(x, b, result);
434
435 // return the result of the solver
436 return result.converged;
437#if HAVE_MPI
438 }
439#endif
440}
441
442template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
444linearSolveBatchwise_(const TracerMatrix& M, std::vector<TracerVector>& x, std::vector<TracerVector>& b)
445{
446 OPM_TIMEBLOCK(tracerSolve);
447 const Scalar tolerance = 1e-2;
448 const int maxIter = 100;
449 const int verbosity = 0;
450
451#if HAVE_MPI
452 if (gridView_.grid().comm().size() > 1)
453 {
454 PropertyTree prm;
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) {
465 x[nrhs] = 0.0;
467 solver->apply(x[nrhs], b[nrhs], result);
468 converged = (converged && result.converged);
469 }
470 return converged;
471 }
472 else
473#endif
474 {
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>;
479
480 if (std::all_of(b.begin(), b.end(),
481 [](const auto& v) { return v.infinity_norm() == 0.0; }))
482 {
483 return true;
484 }
485
486 TracerOperator tracerOperator(M);
487 TracerScalarProduct tracerScalarProduct;
488 TracerPreconditioner tracerPreconditioner(M, 0, 1); // results in ILU0
489
490 TracerSolver solver (tracerOperator, tracerScalarProduct,
491 tracerPreconditioner, tolerance, maxIter,
492 verbosity);
493
494 bool converged = true;
495 for (std::size_t nrhs = 0; nrhs < b.size(); ++nrhs) {
496 x[nrhs] = 0.0;
498 solver.apply(x[nrhs], b[nrhs], result);
499 converged = (converged && result.converged);
500 }
501
502 // return the result of the solver
503 return converged;
504 }
505}
506
507} // namespace Opm
508
509#endif // OPM_GENERIC_TRACER_MODEL_IMPL_HPP
Definition: CollectDataOnIORank.hpp:49
Definition: FlexibleSolver.hpp:45
void setEnableSolTracers(int tracerIdx, bool enableSolTracer)
Definition: GenericTracerModel_impl.hpp:155
int numTracers() const
Return the number of tracers considered by the tracerModel.
Definition: GenericTracerModel_impl.hpp:162
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:227
std::string sname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:176
void setSolTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
Definition: GenericTracerModel_impl.hpp:147
std::string wellfname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:183
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:102
const std::string & name(int tracerIdx) const
Return the tracer name.
Definition: GenericTracerModel_impl.hpp:220
const std::vector< bool > & enableSolTracers() const
Definition: GenericTracerModel_impl.hpp:204
Phase phase(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:197
std::string wellsname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:190
std::string fname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:169
bool linearSolve_(const TracerMatrix &M, TracerVector &x, TracerVector &b)
Definition: GenericTracerModel_impl.hpp:389
Dune::BlockVector< Dune::FieldVector< Scalar, 2 > > TracerVector
Definition: GenericTracerModel.hpp:60
void setFreeTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
Definition: GenericTracerModel_impl.hpp:139
Scalar solTracerConcentration(int tracerIdx, int globalDofIdx) const
Definition: GenericTracerModel_impl.hpp:128
bool linearSolveBatchwise_(const TracerMatrix &M, std::vector< TracerVector > &x, std::vector< TracerVector > &b)
Definition: GenericTracerModel_impl.hpp:444
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 2, 2 > > TracerMatrix
Definition: GenericTracerModel.hpp:59
Scalar freeTracerConcentration(int tracerIdx, int globalDofIdx) const
Return the tracer concentration for tracer index and global DofIdx.
Definition: GenericTracerModel_impl.hpp:117
Scalar currentConcentration_(const Well &eclWell, const std::string &trName, const SummaryState &summaryState) const
Definition: GenericTracerModel_impl.hpp:211
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 TracerSolverSelector< Matrix, Vector >::type > > createParallelFlexibleSolver(const Grid &, const Matrix &, const PropertyTree &)
Definition: GenericTracerModel_impl.hpp:78
Definition: GenericTracerModel_impl.hpp:68
Dune::FlexibleSolver< TracerOperator > type
Definition: GenericTracerModel_impl.hpp:71
Dune::OverlappingSchwarzOperator< M, V, V, Comm > TracerOperator
Definition: GenericTracerModel_impl.hpp:70
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: GenericTracerModel_impl.hpp:69