2021-05-05 05:13:25 -05:00
|
|
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
|
|
|
|
// vi: set et ts=4 sw=4 sts=4:
|
|
|
|
/*
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 2 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
|
|
|
|
Consult the COPYING file in the top-level source directory of this
|
|
|
|
module for the precise wording of the license and the list of
|
|
|
|
copyright holders.
|
|
|
|
*/
|
|
|
|
/**
|
|
|
|
* \file
|
|
|
|
*
|
|
|
|
* \copydoc Opm::EclTracerModel
|
|
|
|
*/
|
|
|
|
#ifndef EWOMS_ECL_GENERIC_TRACER_MODEL_HH
|
|
|
|
#define EWOMS_ECL_GENERIC_TRACER_MODEL_HH
|
|
|
|
|
|
|
|
#include <opm/grid/common/CartesianIndexMapper.hpp>
|
|
|
|
|
|
|
|
#include <opm/models/blackoil/blackoilmodel.hh>
|
|
|
|
#include <opm/common/OpmLog/OpmLog.hpp>
|
|
|
|
|
|
|
|
#include <dune/istl/bcrsmatrix.hh>
|
|
|
|
|
|
|
|
#include <dune/common/version.hh>
|
|
|
|
|
|
|
|
#include <string>
|
|
|
|
#include <vector>
|
|
|
|
#include <iostream>
|
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
class EclipseState;
|
|
|
|
|
|
|
|
template<class Grid, class GridView, class DofMapper, class Stencil, class Scalar>
|
|
|
|
class EclGenericTracerModel {
|
|
|
|
public:
|
2022-06-02 08:12:18 -05:00
|
|
|
using TracerMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;
|
2021-05-05 05:13:25 -05:00
|
|
|
using TracerVector = Dune::BlockVector<Dune::FieldVector<Scalar,1>>;
|
|
|
|
using CartesianIndexMapper = Dune::CartesianIndexMapper<Grid>;
|
2021-09-30 08:13:49 -05:00
|
|
|
static const int dimWorld = Grid::dimensionworld;
|
2021-05-05 05:13:25 -05:00
|
|
|
/*!
|
|
|
|
* \brief Return the number of tracers considered by the tracerModel.
|
|
|
|
*/
|
2021-11-26 10:20:40 -06:00
|
|
|
int numTracers() const;
|
2021-05-05 05:13:25 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Return the tracer name
|
|
|
|
*/
|
2021-11-26 10:21:40 -06:00
|
|
|
const std::string& name(int tracerIdx) const;
|
|
|
|
std::string fname(int tracerIdx) const;
|
|
|
|
|
2021-05-05 05:13:25 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Return the tracer concentration for tracer index and global DofIdx
|
|
|
|
*/
|
|
|
|
Scalar tracerConcentration(int tracerIdx, int globalDofIdx) const;
|
2021-11-23 05:56:56 -06:00
|
|
|
void setTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value);
|
2021-05-05 05:13:25 -05:00
|
|
|
|
2021-06-03 05:58:39 -05:00
|
|
|
/*!
|
|
|
|
* \brief Return well tracer rates
|
|
|
|
*/
|
|
|
|
const std::map<std::pair<std::string, std::string>, double>&
|
|
|
|
getWellTracerRates() const {return wellTracerRate_;}
|
|
|
|
|
2021-05-05 05:13:25 -05:00
|
|
|
protected:
|
|
|
|
EclGenericTracerModel(const GridView& gridView,
|
|
|
|
const EclipseState& eclState,
|
|
|
|
const CartesianIndexMapper& cartMapper,
|
2021-09-30 08:13:49 -05:00
|
|
|
const DofMapper& dofMapper,
|
|
|
|
const std::function<std::array<double,dimWorld>(int)> centroids);
|
2021-05-05 05:13:25 -05:00
|
|
|
|
|
|
|
/*!
|
|
|
|
* \brief Initialize all internal data structures needed by the tracer module
|
|
|
|
*/
|
2021-12-01 03:26:41 -06:00
|
|
|
void doInit(bool rst,
|
2021-05-05 05:13:25 -05:00
|
|
|
size_t numGridDof,
|
|
|
|
size_t gasPhaseIdx,
|
|
|
|
size_t oilPhaseIdx,
|
|
|
|
size_t waterPhaseIdx);
|
|
|
|
|
|
|
|
bool linearSolve_(const TracerMatrix& M, TracerVector& x, TracerVector& b);
|
|
|
|
|
2021-06-03 05:58:39 -05:00
|
|
|
bool linearSolveBatchwise_(const TracerMatrix& M, std::vector<TracerVector>& x, std::vector<TracerVector>& b);
|
|
|
|
|
2021-05-05 05:13:25 -05:00
|
|
|
const GridView& gridView_;
|
|
|
|
const EclipseState& eclState_;
|
|
|
|
const CartesianIndexMapper& cartMapper_;
|
|
|
|
const DofMapper& dofMapper_;
|
|
|
|
|
|
|
|
std::vector<int> tracerPhaseIdx_;
|
|
|
|
std::vector<Dune::BlockVector<Dune::FieldVector<Scalar, 1>>> tracerConcentration_;
|
|
|
|
TracerMatrix *tracerMatrix_;
|
|
|
|
TracerVector tracerResidual_;
|
|
|
|
std::vector<int> cartToGlobal_;
|
|
|
|
std::vector<Dune::BlockVector<Dune::FieldVector<Scalar, 1>>> storageOfTimeIndex1_;
|
2021-06-03 05:58:39 -05:00
|
|
|
|
|
|
|
// <wellName, tracerIdx> -> wellRate
|
|
|
|
std::map<std::pair<std::string, std::string>, double> wellTracerRate_;
|
2021-09-30 08:13:49 -05:00
|
|
|
/// \brief Function returning the cell centers
|
|
|
|
std::function<std::array<double,dimWorld>(int)> centroids_;
|
2021-05-05 05:13:25 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#endif
|