// -*- 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 .
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
#include
#include
#include
#include
#include
#include
#include
namespace Opm {
class EclipseState;
template
class EclGenericTracerModel {
public:
using TracerMatrix = Dune::BCRSMatrix>;
using TracerVector = Dune::BlockVector>;
using CartesianIndexMapper = Dune::CartesianIndexMapper;
static const int dimWorld = Grid::dimensionworld;
/*!
* \brief Return the number of tracers considered by the tracerModel.
*/
int numTracers() const;
/*!
* \brief Return the tracer name
*/
const std::string& name(int tracerIdx) const;
std::string fname(int tracerIdx) const;
/*!
* \brief Return the tracer concentration for tracer index and global DofIdx
*/
Scalar tracerConcentration(int tracerIdx, int globalDofIdx) const;
void setTracerConcentration(int tracerIdx, int globalDofIdx, Scalar value);
/*!
* \brief Return well tracer rates
*/
const std::map, double>&
getWellTracerRates() const {return wellTracerRate_;}
protected:
EclGenericTracerModel(const GridView& gridView,
const EclipseState& eclState,
const CartesianIndexMapper& cartMapper,
const DofMapper& dofMapper,
const std::function(int)> centroids);
/*!
* \brief Initialize all internal data structures needed by the tracer module
*/
void doInit(bool rst,
size_t numGridDof,
size_t gasPhaseIdx,
size_t oilPhaseIdx,
size_t waterPhaseIdx);
bool linearSolve_(const TracerMatrix& M, TracerVector& x, TracerVector& b);
bool linearSolveBatchwise_(const TracerMatrix& M, std::vector& x, std::vector& b);
const GridView& gridView_;
const EclipseState& eclState_;
const CartesianIndexMapper& cartMapper_;
const DofMapper& dofMapper_;
std::vector tracerPhaseIdx_;
std::vector>> tracerConcentration_;
TracerMatrix *tracerMatrix_;
TracerVector tracerResidual_;
std::vector cartToGlobal_;
std::vector>> storageOfTimeIndex1_;
// -> wellRate
std::map, double> wellTracerRate_;
/// \brief Function returning the cell centers
std::function(int)> centroids_;
};
} // namespace Opm
#endif