3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/nonequilibrium/localresidual.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program 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 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program 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 this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
27#define DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
28
29#include <cmath>
32
33namespace Dumux {
34
35// forward declaration
36template<class TypeTag, bool enableChemicalNonEquilibrium>
38
39template <class TypeTag>
41
47template<class TypeTag>
48class NonEquilibriumLocalResidualImplementation<TypeTag, false>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
49{
53 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
54 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
56 using Element = typename GridView::template Codim<0>::Entity;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60
61public:
62 using ParentType::ParentType;
63
73 NumEqVector computeSource(const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const SubControlVolume &scv) const
78 {
79 NumEqVector source(0.0);
80
81 // add contributions from volume flux sources
82 source += problem.source(element, fvGeometry, elemVolVars, scv);
83
84 // add contribution from possible point sources
85 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
86
87 // Call the (kinetic) Energy module, for the source term.
88 // it has to be called from here, because the mass transfered has to be known.
89 if constexpr(ModelTraits::enableThermalNonEquilibrium())
90 {
91 EnergyLocalResidual::computeSourceEnergy(source,
92 element,
93 fvGeometry,
94 elemVolVars,
95 scv);
96 }
97
98 return source;
99 }
100
101};
102
107template<class TypeTag>
108class NonEquilibriumLocalResidualImplementation<TypeTag, true>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
109{
113 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
114 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
115 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
118 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
120 using Element = typename GridView::template Codim<0>::Entity;
121 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
125
127 using Indices = typename ModelTraits::Indices;
128
129 static constexpr int numPhases = ModelTraits::numFluidPhases();
130 static constexpr int numComponents = ModelTraits::numFluidComponents();
131
132 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
133 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
134 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
135 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
136
137 static_assert(numPhases > 1,
138 "chemical non-equlibrium only makes sense for multiple phases");
139 static_assert(numPhases == numComponents,
140 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
141 static_assert(ModelTraits::useMoles(),
142 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
143
144public:
145 using ParentType::ParentType;
153 NumEqVector computeStorage(const Problem& problem,
154 const SubControlVolume& scv,
155 const VolumeVariables& volVars) const
156 {
157 NumEqVector storage(0.0);
158
159 // compute storage term of all components within all phases
160 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
161 {
162 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
163 {
164 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
165 storage[eqIdx] += volVars.porosity()
166 * volVars.saturation(phaseIdx)
167 * volVars.molarDensity(phaseIdx)
168 * volVars.moleFraction(phaseIdx, compIdx);
169 }
171 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
172 }
174 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
175 return storage;
176 }
177
188 NumEqVector computeFlux(const Problem& problem,
189 const Element& element,
190 const FVElementGeometry& fvGeometry,
191 const ElementVolumeVariables& elemVolVars,
192 const SubControlVolumeFace& scvf,
193 const ElementFluxVariablesCache& elemFluxVarsCache) const
194 {
195 FluxVariables fluxVars;
196 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
197 NumEqVector flux(0.0);
198
199 // advective fluxes
200 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
201 {
202 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
203 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
204 {
205 // get equation index
206 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
207
208 // the physical quantities for which we perform upwinding
209 const auto upwindTerm = [phaseIdx, compIdx] (const auto& volVars)
210 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
211
212 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
213
214 // do not add diffusive flux of main component, as that is not done in master as well
215 if (compIdx == phaseIdx)
216 continue;
217 //check for the reference system and adapt units of the diffusive flux accordingly.
218 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
219 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
220 else
221 flux[eqIdx] += diffusiveFluxes[compIdx];
222 }
224 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
225 }
226
228 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
229
230 return flux;
231 }
232
242 NumEqVector computeSource(const Problem& problem,
243 const Element& element,
244 const FVElementGeometry& fvGeometry,
245 const ElementVolumeVariables& elemVolVars,
246 const SubControlVolume &scv) const
247 {
248 NumEqVector source(0.0);
249 // In the case of a kinetic consideration, mass transfer
250 // between phases is realized via source terms there is a
251 // balance equation for each component in each phase
252 const auto& volVars = elemVolVars[scv];
253 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
254
255 //get characteristic length
256 const Scalar characteristicLength = volVars.characteristicLength() ;
257 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
258
259 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
260
261 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
262 {
263 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
264
265 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
266 {
267 if (compIdx <= numPhases)
268 {
269 //we only have to do the source calculation for one component going into the other phase as for each component: n_wn -> n_n = - n_n -> n_wn
270 if (phaseIdx == compIdx)
271 continue;
272
273 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
274
275 //additionally get equilibrium values from volume variables
276 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
277 //get the diffusion coefficient
278 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
279
280 //now compute the flux
281 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
282
283 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
284 componentIntoPhaseMassTransfer[compIdx][compIdx] -= compFluxIntoOtherPhase;
285 }
286 }
287 }
288
289 // Actually add the mass transfer to the sources which might
290 // exist externally
291 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
292 {
293 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
294 {
295 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
296 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
297
298 using std::isfinite;
299 if (!isfinite(source[eqIdx]))
300 DUNE_THROW(NumericalProblem, "Calculated non-finite source");
301 }
302 }
303
304 if constexpr (ModelTraits::enableThermalNonEquilibrium())
305 {
306 // Call the (kinetic) Energy module, for the source term.
307 // it has to be called from here, because the mass transfered has to be known.
308 EnergyLocalResidual::computeSourceEnergy(source,
309 element,
310 fvGeometry,
311 elemVolVars,
312 scv);
313 }
314
315 // add contributions from volume flux sources
316 source += problem.source(element, fvGeometry, elemVolVars, scv);
317
318 // add contribution from possible point sources
319 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
321
322 return source;
323 }
324 };
325
326} // end namespace Dumux
327
328#endif
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
Definition: adapt.hh:29
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:149
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Definition: porousmediumflow/nonequilibrium/localresidual.hh:37
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculates the source term of the equation.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:73
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Calculates the source term of the equation.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:242
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Calculates the storage for all mass balance equations.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:153
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Calculates the flux for all mass balance equations.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:188
Declares all properties used in Dumux.
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...