Skip to content
Snippets Groups Projects
Commit 8e3220f0 authored by Kilian Weishaupt's avatar Kilian Weishaupt Committed by Timo Koch
Browse files

[2p1cni] Use specialized localresidual

parent 0651e76d
No related branches found
No related tags found
Loading
// -*- 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
*
* \brief Element-wise calculation of the residual for problems
* using the two-phase one-component fully implicit models.
*/
#ifndef DUMUX_2P1C_LOCAL_RESIDUAL_HH
#define DUMUX_2P1C_LOCAL_RESIDUAL_HH
#include <dumux/porousmediumflow/immiscible/localresidual.hh>
namespace Dumux
{
/*!
* \ingroup TwoPOneC
* \ingroup ImplicitLocalResidual
* \brief Element-wise calculation of the residual for problems
* using the two-phase one-component fully implicit models.
*/
template<class TypeTag>
class TwoPOneCLocalResidual : public ImmiscibleLocalResidual<TypeTag>
{
using ParentType = ImmiscibleLocalResidual<TypeTag>;
using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
using Problem = typename GET_PROP_TYPE(TypeTag, Problem);
using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
using ResidualVector = typename GET_PROP_TYPE(TypeTag, NumEqVector);
using VolumeVariables = typename GET_PROP_TYPE(TypeTag, VolumeVariables);
using ElementVolumeVariables = typename GET_PROP_TYPE(TypeTag, ElementVolumeVariables);
using FluxVariables = typename GET_PROP_TYPE(TypeTag, FluxVariables);
using ElementFluxVariablesCache = typename GET_PROP_TYPE(TypeTag, ElementFluxVariablesCache);
using SubControlVolume = typename GET_PROP_TYPE(TypeTag, SubControlVolume);
using SubControlVolumeFace = typename GET_PROP_TYPE(TypeTag, SubControlVolumeFace);
using FVElementGeometry = typename GET_PROP_TYPE(TypeTag, FVElementGeometry);
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 Indices = typename GET_PROP_TYPE(TypeTag, Indices);
// first index for the mass balance
enum { conti0EqIdx = Indices::conti0EqIdx };
static const int numPhases = GET_PROP_VALUE(TypeTag, NumPhases);
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 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
*/
ResidualVector computeStorage(const Problem& problem,
const SubControlVolume& scv,
const VolumeVariables& volVars) const
{
ResidualVector storage(0.0);
// compute storage term of all components within all phases
for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
{
storage[conti0EqIdx] +=
volVars.porosity()
* volVars.saturation(phaseIdx) * volVars.density(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 scvf The sub control volume face to compute the flux on
*/
ResidualVector 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);
ResidualVector 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); };
flux[conti0EqIdx] += 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
...@@ -43,6 +43,7 @@ ...@@ -43,6 +43,7 @@
#include "darcyslaw.hh" #include "darcyslaw.hh"
#include "vtkoutputfields.hh" #include "vtkoutputfields.hh"
#include "localresidual.hh"
#include "indices.hh" #include "indices.hh"
#include "volumevariables.hh" #include "volumevariables.hh"
#include "primaryvariableswitch.hh" #include "primaryvariableswitch.hh"
...@@ -132,7 +133,7 @@ namespace Properties ...@@ -132,7 +133,7 @@ namespace Properties
"Only fluid systems with 2 phases are supported by the 2p1cni model!"); "Only fluid systems with 2 phases are supported by the 2p1cni model!");
}; };
SET_TYPE_PROP(TwoPOneCNI, LocalResidual, ImmiscibleLocalResidual<TypeTag>); //! The local residual function SET_TYPE_PROP(TwoPOneCNI, LocalResidual, TwoPOneCLocalResidual<TypeTag>); //! The local residual function
SET_BOOL_PROP(TwoPOneCNI, EnableAdvection, true); //! The one-phase model considers advection SET_BOOL_PROP(TwoPOneCNI, EnableAdvection, true); //! The one-phase model considers advection
SET_BOOL_PROP(TwoPOneCNI, EnableMolecularDiffusion, false); //! The one-phase model has no molecular diffusion SET_BOOL_PROP(TwoPOneCNI, EnableMolecularDiffusion, false); //! The one-phase model has no molecular diffusion
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment