mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-12-02 05:49:09 -06:00
201 lines
7.4 KiB
C++
201 lines
7.4 KiB
C++
// -*- 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::PvsLocalResidual
|
|
*/
|
|
#ifndef EWOMS_PVS_LOCAL_RESIDUAL_HH
|
|
#define EWOMS_PVS_LOCAL_RESIDUAL_HH
|
|
|
|
#include "pvsproperties.hh"
|
|
|
|
#include <opm/models/common/diffusionmodule.hh>
|
|
#include <opm/models/common/energymodule.hh>
|
|
|
|
#include <opm/material/common/Valgrind.hpp>
|
|
|
|
namespace Opm {
|
|
|
|
/*!
|
|
* \ingroup PvsModel
|
|
*
|
|
* \brief Element-wise calculation of the local residual for the
|
|
* compositional multi-phase primary variable switching model.
|
|
*/
|
|
template <class TypeTag>
|
|
class PvsLocalResidual : public GetPropType<TypeTag, Properties::DiscLocalResidual>
|
|
{
|
|
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
|
|
using EqVector = GetPropType<TypeTag, Properties::EqVector>;
|
|
using RateVector = GetPropType<TypeTag, Properties::RateVector>;
|
|
using Indices = GetPropType<TypeTag, Properties::Indices>;
|
|
using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
|
|
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
|
|
|
|
enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
|
|
enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
|
|
enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
|
|
enum { conti0EqIdx = Indices::conti0EqIdx };
|
|
|
|
enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
|
|
using DiffusionModule = Opm::DiffusionModule<TypeTag, enableDiffusion>;
|
|
|
|
enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
|
|
using EnergyModule = Opm::EnergyModule<TypeTag, enableEnergy>;
|
|
|
|
using Toolbox = Opm::MathToolbox<Evaluation>;
|
|
|
|
public:
|
|
/*!
|
|
* \copydoc ImmiscibleLocalResidual::addPhaseStorage
|
|
*/
|
|
template <class LhsEval>
|
|
void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
|
|
const ElementContext& elemCtx,
|
|
unsigned dofIdx,
|
|
unsigned timeIdx,
|
|
unsigned phaseIdx) const
|
|
{
|
|
const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
|
|
const auto& fs = intQuants.fluidState();
|
|
|
|
// compute storage term of all components within all phases
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
|
|
unsigned eqIdx = conti0EqIdx + compIdx;
|
|
storage[eqIdx] +=
|
|
Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx))
|
|
* Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
|
|
* Toolbox::template decay<LhsEval>(intQuants.porosity());
|
|
}
|
|
|
|
EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
|
|
}
|
|
|
|
/*!
|
|
* \copydoc ImmiscibleLocalResidual::computeStorage
|
|
*/
|
|
template <class LhsEval>
|
|
void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
|
|
const ElementContext& elemCtx,
|
|
unsigned dofIdx,
|
|
unsigned timeIdx) const
|
|
{
|
|
storage = 0.0;
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
|
|
addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
|
|
|
|
EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
|
|
}
|
|
|
|
/*!
|
|
* \copydoc ImmiscibleLocalResidual::computeFlux
|
|
*/
|
|
void computeFlux(RateVector& flux,
|
|
const ElementContext& elemCtx,
|
|
unsigned scvfIdx,
|
|
unsigned timeIdx) const
|
|
{
|
|
flux = 0.0;
|
|
addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
|
|
Opm::Valgrind::CheckDefined(flux);
|
|
|
|
addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
|
|
Opm::Valgrind::CheckDefined(flux);
|
|
}
|
|
|
|
/*!
|
|
* \copydoc ImmiscibleLocalResidual::addAdvectiveFlux
|
|
*/
|
|
void addAdvectiveFlux(RateVector& flux,
|
|
const ElementContext& elemCtx,
|
|
unsigned scvfIdx,
|
|
unsigned timeIdx) const
|
|
{
|
|
const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
|
|
|
|
unsigned focusDofIdx = elemCtx.focusDofIndex();
|
|
for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
|
|
// data attached to upstream and the downstream DOFs
|
|
// of the current phase
|
|
unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
|
|
const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
|
|
|
|
// this is a bit hacky because it is specific to the element-centered
|
|
// finite volume scheme. (N.B. that if finite differences are used to
|
|
// linearize the system of equations, it does not matter.)
|
|
if (upIdx == focusDofIdx) {
|
|
Evaluation tmp =
|
|
up.fluidState().molarDensity(phaseIdx)
|
|
* extQuants.volumeFlux(phaseIdx);
|
|
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
|
|
flux[conti0EqIdx + compIdx] +=
|
|
tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
|
|
}
|
|
}
|
|
else {
|
|
Evaluation tmp =
|
|
Toolbox::value(up.fluidState().molarDensity(phaseIdx))
|
|
* extQuants.volumeFlux(phaseIdx);
|
|
|
|
for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
|
|
flux[conti0EqIdx + compIdx] +=
|
|
tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
|
|
}
|
|
}
|
|
}
|
|
|
|
EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
|
|
}
|
|
|
|
/*!
|
|
* \copydoc ImmiscibleLocalResidual::addDiffusiveFlux
|
|
*/
|
|
void addDiffusiveFlux(RateVector& flux,
|
|
const ElementContext& elemCtx,
|
|
unsigned scvfIdx,
|
|
unsigned timeIdx) const
|
|
{
|
|
DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
|
|
EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
|
|
}
|
|
|
|
/*!
|
|
* \copydoc ImmiscibleLocalResidual::computeSource
|
|
*/
|
|
void computeSource(RateVector& source,
|
|
const ElementContext& elemCtx,
|
|
unsigned dofIdx,
|
|
unsigned timeIdx) const
|
|
{
|
|
Opm::Valgrind::SetUndefined(source);
|
|
elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
|
|
Opm::Valgrind::CheckDefined(source);
|
|
}
|
|
};
|
|
|
|
} // namespace Opm
|
|
|
|
#endif
|