GenericTracerModel.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_HPP
29#define OPM_GENERIC_TRACER_MODEL_HPP
30
31#include <dune/istl/bcrsmatrix.hh>
32
33#include <opm/grid/common/CartesianIndexMapper.hpp>
34
35#include <opm/models/blackoil/blackoilmodel.hh>
36
37#include <opm/simulators/linalg/matrixblock.hh>
38
39#include <array>
40#include <cstddef>
41#include <functional>
42#include <map>
43#include <memory>
44#include <string>
45#include <vector>
46
47namespace Opm {
48
49class EclipseState;
50class Well;
51
52template<class Grid, class GridView, class DofMapper, class Stencil, class Scalar>
54public:
55 using TracerMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
56 using TracerVector = Dune::BlockVector<Dune::FieldVector<Scalar,1>>;
58 static constexpr int dimWorld = Grid::dimensionworld;
62 int numTracers() const;
63
67 const std::string& name(int tracerIdx) const;
68 std::string fname(int tracerIdx) const;
69
70
74 Scalar tracerConcentration(int tracerIdx, int globalDofIdx) const;
75 void setTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value);
76
80 const std::map<std::pair<std::string, std::string>, double>&
82
83 template<class Serializer>
84 void serializeOp(Serializer& serializer)
85 {
86 serializer(tracerConcentration_);
87 serializer(wellTracerRate_);
88 }
89
90protected:
91 GenericTracerModel(const GridView& gridView,
92 const EclipseState& eclState,
93 const CartesianIndexMapper& cartMapper,
94 const DofMapper& dofMapper,
95 const std::function<std::array<double,dimWorld>(int)> centroids);
96
100 void doInit(bool rst,
101 std::size_t numGridDof,
102 std::size_t gasPhaseIdx,
103 std::size_t oilPhaseIdx,
104 std::size_t waterPhaseIdx);
105
106 bool linearSolve_(const TracerMatrix& M, TracerVector& x, TracerVector& b);
107
108 bool linearSolveBatchwise_(const TracerMatrix& M, std::vector<TracerVector>& x, std::vector<TracerVector>& b);
109
110 double currentConcentration_(const Well& eclWell, const std::string& name) const;
111
112 const GridView& gridView_;
113 const EclipseState& eclState_;
115 const DofMapper& dofMapper_;
116
117 std::vector<int> tracerPhaseIdx_;
118 std::vector<Dune::BlockVector<Dune::FieldVector<Scalar, 1>>> tracerConcentration_;
119 std::unique_ptr<TracerMatrix> tracerMatrix_;
120 std::vector<Dune::BlockVector<Dune::FieldVector<Scalar, 1>>> storageOfTimeIndex1_;
121
122 // <wellName, tracerIdx> -> wellRate
123 std::map<std::pair<std::string, std::string>, double> wellTracerRate_;
125 std::function<std::array<double,dimWorld>(int)> centroids_;
126};
127
128} // namespace Opm
129
130#endif // OPM_GENERIC_TRACER_MODEL_HPP
Definition: CollectDataOnIORank.hpp:49
Definition: GenericTracerModel.hpp:53
const std::map< std::pair< std::string, std::string >, double > & getWellTracerRates() const
Return well tracer rates.
Definition: GenericTracerModel.hpp:81
Scalar tracerConcentration(int tracerIdx, int globalDofIdx) const
Return the tracer concentration for tracer index and global DofIdx.
Definition: GenericTracerModel_impl.hpp:117
std::vector< Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > > tracerConcentration_
Definition: GenericTracerModel.hpp:118
const GridView & gridView_
Definition: GenericTracerModel.hpp:112
bool linearSolve_(const TracerMatrix &M, TracerVector &x, TracerVector &b)
Definition: GenericTracerModel_impl.hpp:258
std::string fname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:141
static constexpr int dimWorld
Definition: GenericTracerModel.hpp:58
std::unique_ptr< TracerMatrix > tracerMatrix_
Definition: GenericTracerModel.hpp:119
std::vector< int > tracerPhaseIdx_
Definition: GenericTracerModel.hpp:117
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:162
const std::string & name(int tracerIdx) const
Return the tracer name.
Definition: GenericTracerModel_impl.hpp:155
std::function< std::array< double, dimWorld >(int)> centroids_
Function returning the cell centers.
Definition: GenericTracerModel.hpp:125
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > TracerVector
Definition: GenericTracerModel.hpp:56
void setTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
Definition: GenericTracerModel_impl.hpp:127
const DofMapper & dofMapper_
Definition: GenericTracerModel.hpp:115
const CartesianIndexMapper & cartMapper_
Definition: GenericTracerModel.hpp:114
const EclipseState & eclState_
Definition: GenericTracerModel.hpp:113
std::map< std::pair< std::string, std::string >, double > wellTracerRate_
Definition: GenericTracerModel.hpp:123
bool linearSolveBatchwise_(const TracerMatrix &M, std::vector< TracerVector > &x, std::vector< TracerVector > &b)
Definition: GenericTracerModel_impl.hpp:313
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
int numTracers() const
Return the number of tracers considered by the tracerModel.
Definition: GenericTracerModel_impl.hpp:134
Dune::BCRSMatrix< Opm::MatrixBlock< Scalar, 1, 1 > > TracerMatrix
Definition: GenericTracerModel.hpp:55
void serializeOp(Serializer &serializer)
Definition: GenericTracerModel.hpp:84
double currentConcentration_(const Well &eclWell, const std::string &name) const
Definition: GenericTracerModel_impl.hpp:148
std::vector< Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > > storageOfTimeIndex1_
Definition: GenericTracerModel.hpp:120
Dune::CartesianIndexMapper< Grid > CartesianIndexMapper
Definition: GenericTracerModel.hpp:57
Definition: BlackoilPhases.hpp:27