ptflash/flashlocalresidual.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 OPM_PTFLASH_LOCAL_RESIDUAL_HH
29#define OPM_PTFLASH_LOCAL_RESIDUAL_HH
30
31#include <opm/material/common/MathToolbox.hpp>
32#include <opm/material/common/Valgrind.hpp>
33
36
37namespace Opm {
38
45template <class TypeTag>
46class FlashLocalResidual: public GetPropType<TypeTag, Properties::DiscLocalResidual>
47{
48 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
49 using EqVector = GetPropType<TypeTag, Properties::EqVector>;
50 using RateVector = GetPropType<TypeTag, Properties::RateVector>;
51 using Indices = GetPropType<TypeTag, Properties::Indices>;
52 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
53 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
54 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
55
56 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
57 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
58 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
59 enum { water0Idx = Indices::water0Idx };
60 enum { conti0EqIdx = Indices::conti0EqIdx };
61
62 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
63
64 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
66
67 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
69
70 static constexpr bool waterEnabled = Indices::waterEnabled;
71
72 using Toolbox = MathToolbox<Evaluation>;
73
74public:
78 template <class LhsEval>
79 void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
80 const ElementContext& elemCtx,
81 unsigned dofIdx,
82 unsigned timeIdx,
83 unsigned phaseIdx) const
84 {
85 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
86 const auto& fs = intQuants.fluidState();
87
88 // compute water storage term
89 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)) {
90 const unsigned eqIdx = conti0EqIdx + numComponents;
91 storage[eqIdx] =
92 Toolbox::template decay<LhsEval>(fs.density(phaseIdx)) *
93 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
94 Toolbox::template decay<LhsEval>(intQuants.porosity());
95 }
96 else {
97 // compute storage term of all components within oil/gas phases
98 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
99 const unsigned eqIdx = conti0EqIdx + compIdx;
100 storage[eqIdx] +=
101 Toolbox::template decay<LhsEval>(fs.massFraction(phaseIdx, compIdx)) *
102 Toolbox::template decay<LhsEval>(fs.density(phaseIdx)) *
103 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
104 Toolbox::template decay<LhsEval>(intQuants.porosity());
105 }
106 }
107
108 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
109 }
110
114 template <class LhsEval>
115 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
116 const ElementContext& elemCtx,
117 unsigned dofIdx,
118 unsigned timeIdx) const
119 {
120 storage = 0;
121 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
122 addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
123 }
124
125 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
126 }
127
131 void computeFlux(RateVector& flux,
132 const ElementContext& elemCtx,
133 unsigned scvfIdx,
134 unsigned timeIdx) const
135 {
136 flux = 0.0;
137 addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
138 Valgrind::CheckDefined(flux);
139
140 addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
141 Valgrind::CheckDefined(flux);
142 }
143
147 void addAdvectiveFlux(RateVector& flux,
148 const ElementContext& elemCtx,
149 unsigned scvfIdx,
150 unsigned timeIdx) const
151 {
152 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
153
154 const unsigned focusDofIdx = elemCtx.focusDofIndex();
155 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
156 // data attached to upstream and the finite volume of the current phase
157 const auto upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
158 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
159
160 // this is a bit hacky because it is specific to the element-centered
161 // finite volume scheme. (N.B. that if finite differences are used to
162 // linearize the system of equations, it does not matter.)
163 if (upIdx == focusDofIdx) {
164 const Evaluation tmp =
165 up.fluidState().density(phaseIdx) *
166 extQuants.volumeFlux(phaseIdx);
167
168 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)) {
169 const unsigned eqIdx = conti0EqIdx + numComponents;
170 flux[eqIdx] = tmp;
171 }
172 else {
173 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
174 flux[conti0EqIdx + compIdx] +=
175 tmp * up.fluidState().massFraction(phaseIdx, compIdx);
176 }
177 }
178 }
179 else {
180 const Evaluation tmp =
181 Toolbox::value(up.fluidState().density(phaseIdx)) *
182 extQuants.volumeFlux(phaseIdx);
183
184 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)) {
185 const unsigned eqIdx = conti0EqIdx + numComponents;
186 flux[eqIdx] = tmp;
187 }
188 else {
189 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
190 flux[conti0EqIdx + compIdx] +=
191 tmp * Toolbox::value(up.fluidState().massFraction(phaseIdx, compIdx));
192 }
193 }
194 }
195 }
196
197 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
198 }
199
203 void addDiffusiveFlux(RateVector& flux,
204 const ElementContext& elemCtx,
205 unsigned scvfIdx,
206 unsigned timeIdx) const
207 {
208 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
209 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
210 }
211
215 void computeSource(RateVector& source,
216 const ElementContext& elemCtx,
217 unsigned dofIdx,
218 unsigned timeIdx) const
219 {
220 Valgrind::SetUndefined(source);
221 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
222 Valgrind::CheckDefined(source);
223 }
224};
225
226} // namespace Opm
227
228#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
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: ptflash/flashlocalresidual.hh:115
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: flash/flashlocalresidual.hh:174
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: ptflash/flashlocalresidual.hh:79
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: ptflash/flashlocalresidual.hh:215
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: flash/flashlocalresidual.hh:130
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: ptflash/flashlocalresidual.hh:131
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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