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
36
38
39#include <opm/input/eclipse/EclipseState/Phase.hpp>
40
41#include <array>
42#include <cstddef>
43#include <functional>
44#include <map>
45#include <memory>
46#include <string>
47#include <vector>
48
49namespace Opm {
50
51class EclipseState;
52class Well;
53
54template<class Grid, class GridView, class DofMapper, class Stencil, class FluidSystem, class Scalar>
56public:
57 using TracerVectorSingle = Dune::BlockVector<Dune::FieldVector<Scalar, 1>>;
58 using TracerMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 2, 2>>;
59 using TracerVector = Dune::BlockVector<Dune::FieldVector<Scalar, 2>>;
61 static constexpr int dimWorld = Grid::dimensionworld;
65 int numTracers() const;
66
70 const std::string& name(int tracerIdx) const;
71 std::string fname(int tracerIdx) const;
72 std::string sname(int tracerIdx) const;
73 std::string wellfname(int tracerIdx) const;
74 std::string wellsname(int tracerIdx) const;
75
76 Phase phase(int tracerIdx) const;
77 const std::vector<bool>& enableSolTracers() const;
78
82 Scalar freeTracerConcentration(int tracerIdx, int globalDofIdx) const;
83 Scalar solTracerConcentration(int tracerIdx, int globalDofIdx) const;
84 void setFreeTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value);
85 void setSolTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value);
86 void setEnableSolTracers(int tracerIdx, bool enableSolTracer);
87
91 const std::map<std::pair<std::string, std::string>, Scalar>&
93 const std::map<std::pair<std::string, std::string>, Scalar>&
95 const std::map<std::pair<std::string, std::string>, Scalar>&
97 const std::map<std::tuple<std::string, std::string, std::size_t>, Scalar>&
99
100 template<class Serializer>
101 void serializeOp(Serializer& serializer)
102 {
103 serializer(tracerConcentration_);
104 serializer(freeTracerConcentration_);
105 serializer(solTracerConcentration_);
106 serializer(wellTracerRate_);
107 serializer(wellFreeTracerRate_);
108 serializer(wellSolTracerRate_);
109 serializer(mSwTracerRate_);
110 }
111
112protected:
113 GenericTracerModel(const GridView& gridView,
114 const EclipseState& eclState,
115 const CartesianIndexMapper& cartMapper,
116 const DofMapper& dofMapper,
117 const std::function<std::array<double,dimWorld>(int)> centroids);
118
122 void doInit(bool rst,
123 std::size_t numGridDof,
124 std::size_t gasPhaseIdx,
125 std::size_t oilPhaseIdx,
126 std::size_t waterPhaseIdx);
127
128 bool linearSolve_(const TracerMatrix& M, TracerVector& x, TracerVector& b);
129
131 std::vector<TracerVector>& x,
132 std::vector<TracerVector>& b);
133
134 Scalar currentConcentration_(const Well& eclWell, const std::string& name) const;
135
136 const GridView& gridView_;
137 const EclipseState& eclState_;
139 const DofMapper& dofMapper_;
140
141 std::vector<int> tracerPhaseIdx_;
142 std::vector<bool> enableSolTracers_;
143 std::vector<TracerVector> tracerConcentration_;
144 std::unique_ptr<TracerMatrix> tracerMatrix_;
145 std::vector<TracerVectorSingle> freeTracerConcentration_;
146 std::vector<TracerVectorSingle> solTracerConcentration_;
147
148 // <wellName, tracerName> -> wellRate
149 std::map<std::pair<std::string, std::string>, Scalar> wellTracerRate_;
150 std::map<std::pair<std::string, std::string>, Scalar> wellFreeTracerRate_;
151 std::map<std::pair<std::string, std::string>, Scalar> wellSolTracerRate_;
152
153 // <wellName, tracerName, segNum> -> wellRate
154 std::map<std::tuple<std::string, std::string, std::size_t>, Scalar> mSwTracerRate_;
155
157 std::function<std::array<double,dimWorld>(int)> centroids_;
158};
159
160} // namespace Opm
161
162#endif // OPM_GENERIC_TRACER_MODEL_HPP
Definition: CollectDataOnIORank.hpp:49
Definition: GenericTracerModel.hpp:55
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
const std::map< std::pair< std::string, std::string >, Scalar > & getWellSolTracerRates() const
Definition: GenericTracerModel.hpp:96
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
std::vector< TracerVector > tracerConcentration_
Definition: GenericTracerModel.hpp:143
std::vector< int > tracerPhaseIdx_
Definition: GenericTracerModel.hpp:141
std::map< std::tuple< std::string, std::string, std::size_t >, Scalar > mSwTracerRate_
Definition: GenericTracerModel.hpp:154
const DofMapper & dofMapper_
Definition: GenericTracerModel.hpp:139
const GridView & gridView_
Definition: GenericTracerModel.hpp:136
const EclipseState & eclState_
Definition: GenericTracerModel.hpp:137
void setSolTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
Definition: GenericTracerModel_impl.hpp:144
std::string wellfname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:180
std::map< std::pair< std::string, std::string >, Scalar > wellTracerRate_
Definition: GenericTracerModel.hpp:149
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::map< std::pair< std::string, std::string >, Scalar > wellSolTracerRate_
Definition: GenericTracerModel.hpp:151
std::string fname(int tracerIdx) const
Definition: GenericTracerModel_impl.hpp:166
std::vector< TracerVectorSingle > freeTracerConcentration_
Definition: GenericTracerModel.hpp:145
std::function< std::array< double, dimWorld >(int)> centroids_
Function returning the cell centers.
Definition: GenericTracerModel.hpp:157
bool linearSolve_(const TracerMatrix &M, TracerVector &x, TracerVector &b)
Definition: GenericTracerModel_impl.hpp:385
const CartesianIndexMapper & cartMapper_
Definition: GenericTracerModel.hpp:138
Dune::BlockVector< Dune::FieldVector< Scalar, 2 > > TracerVector
Definition: GenericTracerModel.hpp:59
Dune::CartesianIndexMapper< Grid > CartesianIndexMapper
Definition: GenericTracerModel.hpp:60
void setFreeTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value)
Definition: GenericTracerModel_impl.hpp:136
const std::map< std::pair< std::string, std::string >, Scalar > & getWellTracerRates() const
Return well tracer rates.
Definition: GenericTracerModel.hpp:92
std::vector< bool > enableSolTracers_
Definition: GenericTracerModel.hpp:142
Scalar solTracerConcentration(int tracerIdx, int globalDofIdx) const
Definition: GenericTracerModel_impl.hpp:126
const std::map< std::tuple< std::string, std::string, std::size_t >, Scalar > & getMswTracerRates() const
Definition: GenericTracerModel.hpp:98
static constexpr int dimWorld
Definition: GenericTracerModel.hpp:61
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
const std::map< std::pair< std::string, std::string >, Scalar > & getWellFreeTracerRates() const
Definition: GenericTracerModel.hpp:94
void serializeOp(Serializer &serializer)
Definition: GenericTracerModel.hpp:101
std::vector< TracerVectorSingle > solTracerConcentration_
Definition: GenericTracerModel.hpp:146
Scalar freeTracerConcentration(int tracerIdx, int globalDofIdx) const
Return the tracer concentration for tracer index and global DofIdx.
Definition: GenericTracerModel_impl.hpp:116
Dune::BlockVector< Dune::FieldVector< Scalar, 1 > > TracerVectorSingle
Definition: GenericTracerModel.hpp:57
std::map< std::pair< std::string, std::string >, Scalar > wellFreeTracerRate_
Definition: GenericTracerModel.hpp:150
std::unique_ptr< TracerMatrix > tracerMatrix_
Definition: GenericTracerModel.hpp:144
Definition: blackoilboundaryratevector.hh:37