flash/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 EWOMS_FLASH_LOCAL_RESIDUAL_HH
29#define EWOMS_FLASH_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{
54
55 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
56 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
57 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
58 enum { conti0EqIdx = Indices::conti0EqIdx };
59
60 enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
62
63 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
65
66 using Toolbox = Opm::MathToolbox<Evaluation>;
67
68public:
72 template <class LhsEval>
73 void addPhaseStorage(Dune::FieldVector<LhsEval, numEq>& storage,
74 const ElementContext& elemCtx,
75 unsigned dofIdx,
76 unsigned timeIdx,
77 unsigned phaseIdx) const
78 {
79 const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
80 const auto& fs = intQuants.fluidState();
81
82 // compute storage term of all components within all phases
83 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
84 unsigned eqIdx = conti0EqIdx + compIdx;
85 storage[eqIdx] +=
86 Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx)) *
87 Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx)) *
88 Toolbox::template decay<LhsEval>(intQuants.porosity());
89 }
90
91 EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
92 }
93
97 template <class LhsEval>
98 void computeStorage(Dune::FieldVector<LhsEval, numEq>& storage,
99 const ElementContext& elemCtx,
100 unsigned dofIdx,
101 unsigned timeIdx) const
102 {
103 storage = 0;
104 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
105 addPhaseStorage(storage, elemCtx, dofIdx, timeIdx, phaseIdx);
106 }
107
108 EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
109 }
110
114 void computeFlux(RateVector& flux,
115 const ElementContext& elemCtx,
116 unsigned scvfIdx,
117 unsigned timeIdx) const
118 {
119 flux = 0.0;
120 addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
121 Valgrind::CheckDefined(flux);
122
123 addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
124 Valgrind::CheckDefined(flux);
125 }
126
130 void addAdvectiveFlux(RateVector& flux,
131 const ElementContext& elemCtx,
132 unsigned scvfIdx,
133 unsigned timeIdx) const
134 {
135 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
136
137 const unsigned focusDofIdx = elemCtx.focusDofIndex();
138 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
139 // data attached to upstream and the finite volume of the current phase
140 const unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
141 const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
142
143 // this is a bit hacky because it is specific to the element-centered
144 // finite volume scheme. (N.B. that if finite differences are used to
145 // linearize the system of equations, it does not matter.)
146 if (upIdx == focusDofIdx) {
147 Evaluation tmp =
148 up.fluidState().molarDensity(phaseIdx) *
149 extQuants.volumeFlux(phaseIdx);
150
151 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
152 flux[conti0EqIdx + compIdx] +=
153 tmp * up.fluidState().moleFraction(phaseIdx, compIdx);
154 }
155 }
156 else {
157 Evaluation tmp =
158 Toolbox::value(up.fluidState().molarDensity(phaseIdx)) *
159 extQuants.volumeFlux(phaseIdx);
160
161 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
162 flux[conti0EqIdx + compIdx] +=
163 tmp * Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
164 }
165 }
166 }
167
168 EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
169 }
170
174 void addDiffusiveFlux(RateVector& flux,
175 const ElementContext& elemCtx,
176 unsigned scvfIdx,
177 unsigned timeIdx) const
178 {
179 DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
180 EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
181 }
182
186 void computeSource(RateVector& source,
187 const ElementContext& elemCtx,
188 unsigned dofIdx,
189 unsigned timeIdx) const
190 {
191 Valgrind::SetUndefined(source);
192 elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
193 Valgrind::CheckDefined(source);
194 }
195};
196
197} // namespace Opm
198
199#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
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Definition: flash/flashlocalresidual.hh:47
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: flash/flashlocalresidual.hh:98
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: flash/flashlocalresidual.hh:73
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: flash/flashlocalresidual.hh:186
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: flash/flashlocalresidual.hh:114
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