My Project
Loading...
Searching...
No Matches
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
33
34#include <opm/material/common/Valgrind.hpp>
35
36namespace Opm {
43template <class TypeTag>
44class FlashLocalResidual: public GetPropType<TypeTag, Properties::DiscLocalResidual>
45{
53
56 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
57 enum { water0Idx = Indices::water0Idx };
58 enum { conti0EqIdx = Indices::conti0EqIdx };
59
60 enum { waterPhaseIdx = FluidSystem::waterPhaseIdx };
61
62 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
64
67
68 static const bool waterEnabled = Indices::waterEnabled;
69
70 using Toolbox = Opm::MathToolbox<Evaluation>;
71
72public:
76 template <class LhsEval>
77 void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
78 const ElementContext& elemCtx,
79 unsigned dofIdx,
80 unsigned timeIdx,
81 unsigned phaseIdx) const
82 {
83 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
84 const auto& fs = intQuants.fluidState();
85
86 // compute water storage term
87 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)) {
88 unsigned eqIdx = conti0EqIdx + numComponents;
89 storage[eqIdx] =
90 Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
91 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
92 * Toolbox::template decay<LhsEval>(intQuants.porosity());
93 }
94 else {
95 // compute storage term of all components within oil/gas phases
96 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
97 unsigned eqIdx = conti0EqIdx + compIdx;
98 storage[eqIdx] +=
99 Toolbox::template decay<LhsEval>(fs.massFraction(phaseIdx, compIdx))
100 * Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
101 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
102 * Toolbox::template decay<LhsEval>(intQuants.porosity());
103 }
104 }
105
106 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
107 }
108
112 template <class LhsEval>
113 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
114 const ElementContext& elemCtx,
115 unsigned dofIdx,
116 unsigned timeIdx) const
117 {
118 storage = 0;
119 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
120 addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
121
122 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
123 }
124
128 void computeFlux(RateVector& flux,
129 const ElementContext& elemCtx,
130 unsigned scvfIdx,
131 unsigned timeIdx) const
132 {
133 flux = 0.0;
135 Opm::Valgrind::CheckDefined(flux);
136
138 Opm::Valgrind::CheckDefined(flux);
139 }
140
144 void addAdvectiveFlux(RateVector& flux,
145 const ElementContext& elemCtx,
146 unsigned scvfIdx,
147 unsigned timeIdx) const
148 {
149 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
150
151 unsigned focusDofIdx = elemCtx.focusDofIndex();
152 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
153 // data attached to upstream and the finite volume of the current phase
154 auto upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
155 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
156
157 // this is a bit hacky because it is specific to the element-centered
158 // finite volume scheme. (N.B. that if finite differences are used to
159 // linearize the system of equations, it does not matter.)
160 if (upIdx == focusDofIdx) {
161 Evaluation tmp =
162 up.fluidState().density(phaseIdx)
163 * extQuants.volumeFlux(phaseIdx);
164
165 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)) {
166 unsigned eqIdx = conti0EqIdx + numComponents;
167 flux[eqIdx] = tmp;
168 }
169 else {
170 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
171 flux[conti0EqIdx + compIdx] +=
172 tmp*up.fluidState().massFraction(phaseIdx, compIdx);
173 }
174 }
175 }
176 else {
177 Evaluation tmp =
178 Toolbox::value(up.fluidState().density(phaseIdx))
179 * extQuants.volumeFlux(phaseIdx);
180
181 if (waterEnabled && phaseIdx == static_cast<unsigned int>(waterPhaseIdx)) {
182 unsigned eqIdx = conti0EqIdx + numComponents;
183 flux[eqIdx] = tmp;
184 }
185 else {
186 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
187 flux[conti0EqIdx + compIdx] +=
188 tmp*Toolbox::value(up.fluidState().massFraction(phaseIdx, compIdx));
189 }
190 }
191 }
192 }
193
194 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
195 }
196
200 void addDiffusiveFlux(RateVector& flux,
201 const ElementContext& elemCtx,
202 unsigned scvfIdx,
203 unsigned timeIdx) const
204 {
205 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
206 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
207 }
208
212 void computeSource(RateVector& source,
213 const ElementContext& elemCtx,
214 unsigned dofIdx,
215 unsigned timeIdx) const
216 {
217 Opm::Valgrind::SetUndefined(source);
218 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
219 Opm::Valgrind::CheckDefined(source);
220 }
221};
222
223} // namespace Opm
224
225#endif
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition diffusionmodule.hh:48
Provides the auxiliary methods required for consideration of the energy equation.
Definition energymodule.hh:50
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Evaluate the amount all conservation quantities (e.g.
Definition flashlocalresidual.hh:113
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition flashlocalresidual.hh:172
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.
Definition flashlocalresidual.hh:77
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition flashlocalresidual.hh:212
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition flashlocalresidual.hh:128
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 flashlocalresidual.hh:128
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
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:235