ncplocalresidual.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_NCP_LOCAL_RESIDUAL_HH
29#define EWOMS_NCP_LOCAL_RESIDUAL_HH
30
31#include <dune/istl/bvector.hh>
32
33#include <opm/material/common/MathToolbox.hpp>
34#include <opm/material/common/Valgrind.hpp>
35
40
41#include <type_traits>
42
43namespace Opm {
44
51template <class TypeTag>
52class NcpLocalResidual : public GetPropType<TypeTag, Properties::DiscLocalResidual>
53{
62
63 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
64 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
65 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
66 enum { ncp0EqIdx = Indices::ncp0EqIdx };
67 enum { conti0EqIdx = Indices::conti0EqIdx };
68
69 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
71
72 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
74
75 using EvalEqVector = Dune::FieldVector<Evaluation, numEq>;
76 using ElemEvalEqVector = Dune::BlockVector<EvalEqVector>;
77 using Toolbox = MathToolbox<Evaluation>;
78
79public:
83 template <class LhsEval>
84 void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
85 const ElementContext& elemCtx,
86 unsigned dofIdx,
87 unsigned timeIdx,
88 unsigned phaseIdx) const
89 {
90 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
91 const auto& fluidState = intQuants.fluidState();
92
93 // compute storage term of all components within all phases
94 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
95 const unsigned eqIdx = conti0EqIdx + compIdx;
96 storage[eqIdx] +=
97 Toolbox::template decay<LhsEval>(fluidState.molarity(phaseIdx, compIdx)) *
98 Toolbox::template decay<LhsEval>(fluidState.saturation(phaseIdx)) *
99 Toolbox::template decay<LhsEval>(intQuants.porosity());
100 }
101
102 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
103 }
104
108 template <class LhsEval>
109 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
110 const ElementContext& elemCtx,
111 unsigned dofIdx,
112 unsigned timeIdx) const
113 {
114 storage = 0;
115 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
116 addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
117 }
118
119 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
120 }
121
125 void computeFlux(RateVector& flux,
126 const ElementContext& elemCtx,
127 unsigned scvfIdx,
128 unsigned timeIdx) const
129 {
130 flux = 0.0;
131 addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
132 Valgrind::CheckDefined(flux);
133
134 addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
135 Valgrind::CheckDefined(flux);
136 }
137
141 void addAdvectiveFlux(RateVector& flux,
142 const ElementContext& elemCtx,
143 unsigned scvfIdx,
144 unsigned timeIdx) const
145 {
146 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
147
148 unsigned focusDofIdx = elemCtx.focusDofIndex();
149 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
150 // data attached to upstream and the downstream DOFs
151 // of the current phase
152 const unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
153 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
154
155 // this is a bit hacky because it is specific to the element-centered
156 // finite volume scheme. (N.B. that if finite differences are used to
157 // linearize the system of equations, it does not matter.)
158 if (upIdx == focusDofIdx) {
159 const Evaluation tmp =
160 up.fluidState().molarDensity(phaseIdx) *
161 extQuants.volumeFlux(phaseIdx);
162
163 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
164 flux[conti0EqIdx + compIdx] +=
165 tmp * up.fluidState().moleFraction(phaseIdx, compIdx);
166 }
167 }
168 else {
169 const Evaluation tmp =
170 Toolbox::value(up.fluidState().molarDensity(phaseIdx)) *
171 extQuants.volumeFlux(phaseIdx);
172
173 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
174 flux[conti0EqIdx + compIdx] +=
175 tmp * Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
176 }
177 }
178 }
179
180 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
181 }
182
186 void addDiffusiveFlux(RateVector& flux,
187 const ElementContext& elemCtx,
188 unsigned scvfIdx,
189 unsigned timeIdx) const
190 {
191 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
192 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
193 }
194
201 void computeSource(RateVector& source,
202 const ElementContext& elemCtx,
203 unsigned dofIdx,
204 unsigned timeIdx) const
205 {
206 Valgrind::SetUndefined(source);
207 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
208 Valgrind::CheckDefined(source);
209
210 // evaluate the NCPs (i.e., the "phase presence" equations)
211 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
212 source[ncp0EqIdx + phaseIdx] = phaseNcp(elemCtx, dofIdx, timeIdx, phaseIdx);
213 }
214 }
215
219 template <class LhsEval = Evaluation>
220 LhsEval phaseNcp(const ElementContext& elemCtx,
221 unsigned dofIdx,
222 unsigned timeIdx,
223 unsigned phaseIdx) const
224 {
225 const auto& fluidState = elemCtx.intensiveQuantities(dofIdx, timeIdx).fluidState();
226 using FluidState = std::remove_const_t<std::remove_reference_t<decltype(fluidState)>>;
227
228 using LhsToolbox = MathToolbox<LhsEval>;
229
230 const LhsEval& a = phaseNotPresentIneq_<FluidState, LhsEval>(fluidState, phaseIdx);
231 const LhsEval& b = phasePresentIneq_<FluidState, LhsEval>(fluidState, phaseIdx);
232 return LhsToolbox::min(a, b);
233 }
234
235private:
240 template <class FluidState, class LhsEval>
241 LhsEval phasePresentIneq_(const FluidState& fluidState, unsigned phaseIdx) const
242 {
243 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
244 return FsToolbox::template decay<LhsEval>(fluidState.saturation(phaseIdx));
245 }
246
251 template <class FluidState, class LhsEval>
252 LhsEval phaseNotPresentIneq_(const FluidState& fluidState, unsigned phaseIdx) const
253 {
254 using FsToolbox = MathToolbox<typename FluidState::Scalar>;
255
256 // difference of sum of mole fractions in the phase from 100%
257 LhsEval a = 1.0;
258 for (unsigned i = 0; i < numComponents; ++i) {
259 a -= FsToolbox::template decay<LhsEval>(fluidState.moleFraction(phaseIdx, i));
260 }
261 return a;
262 }
263};
264
265} // namespace Opm
266
267#endif
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:51
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:54
Details needed to calculate the local residual in the compositional multi-phase NCP-model .
Definition: ncplocalresidual.hh:53
LhsEval phaseNcp(const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Returns the value of the NCP-function for a phase.
Definition: ncplocalresidual.hh:220
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: ncplocalresidual.hh:141
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: ncplocalresidual.hh:186
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: ncplocalresidual.hh:125
void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: ncplocalresidual.hh:84
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: ncplocalresidual.hh:109
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: ncplocalresidual.hh:201
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:39
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:233