-
Bernd Flemisch authoredBernd Flemisch authored
mylocalresidual.hh 6.51 KiB
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program 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. *
* *
* This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup PorousmediumImmiscible
* \brief Element-wise calculation of the residual for problems
* using the n-phase immiscible fully implicit models.
*/
#ifndef DUMUX_MY_LOCAL_RESIDUAL_HH
#define DUMUX_MY_LOCAL_RESIDUAL_HH
#include <dumux/common/properties.hh>
namespace Dumux
{
/*!
* \ingroup PorousmediumImmiscible
* \brief Element-wise calculation of the residual for problems
* using the n-phase immiscible fully implicit models.
*/
template<class TypeTag>
class MyLocalResidual : public GET_PROP_TYPE(TypeTag, BaseLocalResidual)
{
using ParentType = typename GET_PROP_TYPE(TypeTag, BaseLocalResidual);
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
using NumEqVector = typename GET_PROP_TYPE(TypeTag, NumEqVector);
using VolumeVariables = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
using ElementVolumeVariables = typename GET_PROP_TYPE(TypeTag, GridVolumeVariables)::LocalView;
using FluxVariables = typename GET_PROP_TYPE(TypeTag, FluxVariables);
using ElementFluxVariablesCache = typename GET_PROP_TYPE(TypeTag, GridFluxVariablesCache)::LocalView;
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVGridGeometry)::LocalView;
using SubControlVolume = typename FVElementGeometry::SubControlVolume;
using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
using Element = typename GridView::template Codim<0>::Entity;
using EnergyLocalResidual = typename GET_PROP_TYPE(TypeTag, EnergyLocalResidual);
using ModelTraits = typename GET_PROP_TYPE(TypeTag, ModelTraits);
static constexpr int numPhases = ModelTraits::numPhases();
static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx; //!< first index for the mass balance
public:
using ParentType::ParentType;
/*!
* \brief Evaluate the rate of change of all conservation
* quantites (e.g. phase mass) within a sub-control
* volume of a finite volume element for the immiscible models.
* \param problem The problem
* \param scv The sub control volume
* \param volVars The current or previous volVars
* \note This function should not include the source and sink terms.
* \note The volVars can be different to allow computing
* the implicit euler time derivative here
*/
// TODO: eliminate density from the storage term
NumEqVector computeStorage(const Problem& problem,
const SubControlVolume& scv,
const VolumeVariables& volVars) const
{
// partial time derivative of the phase mass
NumEqVector storage;
for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
{
auto eqIdx = conti0EqIdx + phaseIdx;
storage[eqIdx] = volVars.porosity()
* volVars.density(phaseIdx)
* volVars.saturation(phaseIdx);
//! The energy storage in the fluid phase with index phaseIdx
EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
}
//! The energy storage in the solid matrix
EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
return storage;
}
/*!
* \brief Evaluate the mass flux over a face of a sub control volume
*
* \param problem The problem
* \param element The element
* \param fvGeometry The finite volume geometry context
* \param elemVolVars The volume variables for all flux stencil elements
* \param scvf The sub control volume face to compute the flux on
* \param elemFluxVarsCache The cache related to flux compuation
*/
// TODO: eliminate the density from the flux term
NumEqVector computeFlux(const Problem& problem,
const Element& element,
const FVElementGeometry& fvGeometry,
const ElementVolumeVariables& elemVolVars,
const SubControlVolumeFace& scvf,
const ElementFluxVariablesCache& elemFluxVarsCache) const
{
FluxVariables fluxVars;
fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
NumEqVector flux;
for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
{
// the physical quantities for which we perform upwinding
auto upwindTerm = [phaseIdx](const auto& volVars)
{ return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); };
auto eqIdx = conti0EqIdx + phaseIdx;
flux[eqIdx] = fluxVars.advectiveFlux(phaseIdx, upwindTerm);
//! Add advective phase energy fluxes. For isothermal model the contribution is zero.
EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
}
//! Add diffusive energy fluxes. For isothermal model the contribution is zero.
EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
return flux;
}
};
} // end namespace Dumux
#endif