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
33
34#include <opm/material/common/Valgrind.hpp>
35
36namespace Opm {
43template <class TypeTag>
44class FlashLocalResidual: public GetPropType<TypeTag, Properties::DiscLocalResidual>
45{
46 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
47 using EqVector = GetPropType<TypeTag, Properties::EqVector>;
48 using RateVector = GetPropType<TypeTag, Properties::RateVector>;
49 using Indices = GetPropType<TypeTag, Properties::Indices>;
50 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
51 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
52
53 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
54 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
55 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
56 enum { conti0EqIdx = Indices::conti0EqIdx };
57
58 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
60
61 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
63
64 using Toolbox = Opm::MathToolbox<Evaluation>;
65
66public:
70 template <class LhsEval>
71 void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
72 const ElementContext& elemCtx,
73 unsigned dofIdx,
74 unsigned timeIdx,
75 unsigned phaseIdx) const
76 {
77 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
78 const auto& fs = intQuants.fluidState();
79
80 // compute storage term of all components within all phases
81 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
82 unsigned eqIdx = conti0EqIdx + compIdx;
83 storage[eqIdx] +=
84 Toolbox::template decay<LhsEval>(fs.massFraction(phaseIdx, compIdx))
85 * Toolbox::template decay<LhsEval>(fs.density(phaseIdx))
86 * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
87 * Toolbox::template decay<LhsEval>(intQuants.porosity());
88 }
89
90 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
91 }
92
96 template <class LhsEval>
97 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
98 const ElementContext& elemCtx,
99 unsigned dofIdx,
100 unsigned timeIdx) const
101 {
102 storage = 0;
103 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
104 addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
105
106 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
107 }
108
112 void computeFlux(RateVector& flux,
113 const ElementContext& elemCtx,
114 unsigned scvfIdx,
115 unsigned timeIdx) const
116 {
117 flux = 0.0;
118 addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
119 Opm::Valgrind::CheckDefined(flux);
120
121 addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
122 Opm::Valgrind::CheckDefined(flux);
123 }
124
128 void addAdvectiveFlux(RateVector& flux,
129 const ElementContext& elemCtx,
130 unsigned scvfIdx,
131 unsigned timeIdx) const
132 {
133 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
134
135 unsigned focusDofIdx = elemCtx.focusDofIndex();
136 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
137 // data attached to upstream and the finite volume of the current phase
138 auto upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
139 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
140
141 // this is a bit hacky because it is specific to the element-centered
142 // finite volume scheme. (N.B. that if finite differences are used to
143 // linearize the system of equations, it does not matter.)
144 if (upIdx == focusDofIdx) {
145 Evaluation tmp =
146 up.fluidState().density(phaseIdx)
147 * extQuants.volumeFlux(phaseIdx);
148
149 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
150 flux[conti0EqIdx + compIdx] +=
151 tmp*up.fluidState().massFraction(phaseIdx, compIdx);
152 }
153 }
154 else {
155 Evaluation tmp =
156 Toolbox::value(up.fluidState().density(phaseIdx))
157 * extQuants.volumeFlux(phaseIdx);
158
159 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
160 flux[conti0EqIdx + compIdx] +=
161 tmp*Toolbox::value(up.fluidState().massFraction(phaseIdx, compIdx));
162 }
163 }
164 }
165
166 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
167 }
168
172 void addDiffusiveFlux(RateVector& flux,
173 const ElementContext& elemCtx,
174 unsigned scvfIdx,
175 unsigned timeIdx) const
176 {
177 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
178 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
179 }
180
184 void computeSource(RateVector& source,
185 const ElementContext& elemCtx,
186 unsigned dofIdx,
187 unsigned timeIdx) const
188 {
189 Opm::Valgrind::SetUndefined(source);
190 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
191 Opm::Valgrind::CheckDefined(source);
192 }
193};
194
195} // namespace Opm
196
197#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. phase mass) within a finite sub-control volume.
Definition: ptflash/flashlocalresidual.hh:97
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: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. phase mass) within a single fluid phase.
Definition: ptflash/flashlocalresidual.hh:71
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: ptflash/flashlocalresidual.hh:184
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: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: ptflash/flashlocalresidual.hh:112
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
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:235