richardslocalresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_RICHARDS_LOCAL_RESIDUAL_HH
29#define EWOMS_RICHARDS_LOCAL_RESIDUAL_HH
30
32
34
35namespace Opm {
36
41template <class TypeTag>
42class RichardsLocalResidual : public GetPropType<TypeTag, Properties::DiscLocalResidual>
43{
50
51 enum { contiEqIdx = Indices::contiEqIdx };
52 enum { liquidPhaseIdx = getPropValue<TypeTag, Properties::LiquidPhaseIndex>() };
53 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
54 using Toolbox = Opm::MathToolbox<Evaluation>;
55
56public:
60 template <class LhsEval>
61 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
62 const ElementContext& elemCtx,
63 unsigned dofIdx,
64 unsigned timeIdx) const
65 {
66 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
67
68 // partial time derivative of the wetting phase mass
69 storage[contiEqIdx] =
70 Toolbox::template decay<LhsEval>(intQuants.fluidState().density(liquidPhaseIdx))
71 *Toolbox::template decay<LhsEval>(intQuants.fluidState().saturation(liquidPhaseIdx))
72 *Toolbox::template decay<LhsEval>(intQuants.porosity());
73 }
74
78 void computeFlux(RateVector& flux,
79 const ElementContext& elemCtx,
80 unsigned scvfIdx,
81 unsigned timeIdx) const
82 {
83 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
84
85 unsigned focusDofIdx = elemCtx.focusDofIndex();
86 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(liquidPhaseIdx));
87
88 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
89
90 // compute advective mass flux of the liquid phase. This is slightly hacky
91 // because it is specific to the element-centered finite volume method.
92 const Evaluation& rho = up.fluidState().density(liquidPhaseIdx);
93 if (focusDofIdx == upIdx)
94 flux[contiEqIdx] = extQuants.volumeFlux(liquidPhaseIdx)*rho;
95 else
96 flux[contiEqIdx] = extQuants.volumeFlux(liquidPhaseIdx)*Toolbox::value(rho);
97 }
98
102 void computeSource(RateVector& source,
103 const ElementContext& elemCtx,
104 unsigned dofIdx,
105 unsigned timeIdx) const
106 { elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx); }
107};
108
109} // namespace Opm
110
111#endif
Element-wise calculation of the residual for the Richards model.
Definition: richardslocalresidual.hh:43
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: richardslocalresidual.hh:78
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume.
Definition: richardslocalresidual.hh:61
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: richardslocalresidual.hh:102
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:242