pvslocalresidual.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_PVS_LOCAL_RESIDUAL_HH
29#define EWOMS_PVS_LOCAL_RESIDUAL_HH
30
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
33
38
39namespace Opm {
40
47template <class TypeTag>
48class PvsLocalResidual : public GetPropType<TypeTag, Properties::DiscLocalResidual>
49{
56
57 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
58 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
59 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
61
62 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
64
65 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
67
68 using Toolbox = MathToolbox<Evaluation>;
69
70public:
74 template <class LhsEval>
75 void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
76 const ElementContext& elemCtx,
77 unsigned dofIdx,
78 unsigned timeIdx,
79 unsigned phaseIdx) const
80 {
81 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
82 const auto& fs = intQuants.fluidState();
83
84 // compute storage term of all components within all phases
85 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
86 const unsigned eqIdx = conti0EqIdx + compIdx;
87 storage[eqIdx] +=
88 Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx)) *
89 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
90 Toolbox::template decay<LhsEval>(intQuants.porosity());
91 }
92
93 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
94 }
95
99 template <class LhsEval>
100 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
101 const ElementContext& elemCtx,
102 unsigned dofIdx,
103 unsigned timeIdx) const
104 {
105 storage = 0.0;
106 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
107 addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
108 }
109
110 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
111 }
112
116 void computeFlux(RateVector& flux,
117 const ElementContext& elemCtx,
118 unsigned scvfIdx,
119 unsigned timeIdx) const
120 {
121 flux = 0.0;
122 addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
123 Valgrind::CheckDefined(flux);
124
125 addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
126 Valgrind::CheckDefined(flux);
127 }
128
132 void addAdvectiveFlux(RateVector& flux,
133 const ElementContext& elemCtx,
134 unsigned scvfIdx,
135 unsigned timeIdx) const
136 {
137 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
138
139 const unsigned focusDofIdx = elemCtx.focusDofIndex();
140 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
141 // data attached to upstream and the downstream DOFs
142 // of the current phase
143 const unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
144 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
145
146 // this is a bit hacky because it is specific to the element-centered
147 // finite volume scheme. (N.B. that if finite differences are used to
148 // linearize the system of equations, it does not matter.)
149 if (upIdx == focusDofIdx) {
150 const Evaluation tmp =
151 up.fluidState().molarDensity(phaseIdx) *
152 extQuants.volumeFlux(phaseIdx);
153
154 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
155 flux[conti0EqIdx + compIdx] +=
156 tmp * up.fluidState().moleFraction(phaseIdx, compIdx);
157 }
158 }
159 else {
160 const Evaluation tmp =
161 Toolbox::value(up.fluidState().molarDensity(phaseIdx)) *
162 extQuants.volumeFlux(phaseIdx);
163
164 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
165 flux[conti0EqIdx + compIdx] +=
166 tmp * Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
167 }
168 }
169 }
170
171 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
172 }
173
177 void addDiffusiveFlux(RateVector& flux,
178 const ElementContext& elemCtx,
179 unsigned scvfIdx,
180 unsigned timeIdx) const
181 {
182 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
183 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
184 }
185
189 void computeSource(RateVector& source,
190 const ElementContext& elemCtx,
191 unsigned dofIdx,
192 unsigned timeIdx) const
193 {
194 Valgrind::SetUndefined(source);
195 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
196 Valgrind::CheckDefined(source);
197 }
198};
199
200} // namespace Opm
201
202#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
Element-wise calculation of the local residual for the compositional multi-phase primary variable swi...
Definition: pvslocalresidual.hh:49
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: pvslocalresidual.hh:75
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: pvslocalresidual.hh:100
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: pvslocalresidual.hh:189
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: pvslocalresidual.hh:177
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: pvslocalresidual.hh:132
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: pvslocalresidual.hh:116
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