version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
15#define DUMUX_NONEQUILIBRIUM_LOCAL_RESIDUAL_HH
16
17#include <cmath>
21
22namespace Dumux {
23
24// forward declaration
25template<class TypeTag, bool enableChemicalNonEquilibrium>
27
28template <class TypeTag>
30
36template<class TypeTag>
37class NonEquilibriumLocalResidualImplementation<TypeTag, false>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
38{
42 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
43 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
45 using Element = typename GridView::template Codim<0>::Entity;
46 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49
50public:
51 using ParentType::ParentType;
52
62 NumEqVector computeSource(const Problem& problem,
63 const Element& element,
64 const FVElementGeometry& fvGeometry,
65 const ElementVolumeVariables& elemVolVars,
66 const SubControlVolume &scv) const
67 {
68 NumEqVector source(0.0);
69
70 // add contributions from volume flux sources
71 source += problem.source(element, fvGeometry, elemVolVars, scv);
72
73 // add contribution from possible point sources
74 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
75
76 // Call the (kinetic) Energy module, for the source term.
77 // it has to be called from here, because the mass transferred has to be known.
78 if constexpr(ModelTraits::enableThermalNonEquilibrium())
79 {
80 EnergyLocalResidual::computeSourceEnergy(source,
81 element,
82 fvGeometry,
83 elemVolVars,
84 scv);
85 }
86
87 return source;
88 }
89
90};
91
96template<class TypeTag>
97class NonEquilibriumLocalResidualImplementation<TypeTag, true>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
98{
102 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
103 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
104 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
107 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
109 using Element = typename GridView::template Codim<0>::Entity;
110 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
114
116 using Indices = typename ModelTraits::Indices;
117
118 static constexpr int numPhases = ModelTraits::numFluidPhases();
119 static constexpr int numComponents = ModelTraits::numFluidComponents();
120
121 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
122 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
123 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
124 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
125
126 static_assert(numPhases > 1,
127 "chemical non-equlibrium only makes sense for multiple phases");
128 static_assert(numPhases == numComponents,
129 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
130 static_assert(ModelTraits::useMoles(),
131 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
132
133public:
134 using ParentType::ParentType;
142 NumEqVector computeStorage(const Problem& problem,
143 const SubControlVolume& scv,
144 const VolumeVariables& volVars) const
145 {
146 NumEqVector storage(0.0);
147
148 // compute storage term of all components within all phases
149 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
150 {
151 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
152 {
153 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
154 storage[eqIdx] += volVars.porosity()
155 * volVars.saturation(phaseIdx)
156 * volVars.molarDensity(phaseIdx)
157 * volVars.moleFraction(phaseIdx, compIdx);
158 }
160 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
161 }
163 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
164 return storage;
165 }
166
177 NumEqVector computeFlux(const Problem& problem,
178 const Element& element,
179 const FVElementGeometry& fvGeometry,
180 const ElementVolumeVariables& elemVolVars,
181 const SubControlVolumeFace& scvf,
182 const ElementFluxVariablesCache& elemFluxVarsCache) const
183 {
184 FluxVariables fluxVars;
185 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
186 NumEqVector flux(0.0);
187
188 // advective fluxes
189 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
190 {
191 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
192 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
193 {
194 // get equation index
195 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
196
197 // the physical quantities for which we perform upwinding
198 const auto upwindTerm = [phaseIdx, compIdx] (const auto& volVars)
199 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
200
201 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
202
203 // do not add diffusive flux of main component, as that is not done in master as well
204 if (compIdx == phaseIdx)
205 continue;
206 //check for the reference system and adapt units of the diffusive flux accordingly.
207 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
208 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
209 else
210 flux[eqIdx] += diffusiveFluxes[compIdx];
211 }
213 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
214 }
215
217 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
218
219 return flux;
220 }
221
231 NumEqVector computeSource(const Problem& problem,
232 const Element& element,
233 const FVElementGeometry& fvGeometry,
234 const ElementVolumeVariables& elemVolVars,
235 const SubControlVolume &scv) const
236 {
237 NumEqVector source(0.0);
238 // In the case of a kinetic consideration, mass transfer
239 // between phases is realized via source terms there is a
240 // balance equation for each component in each phase
241 const auto& volVars = elemVolVars[scv];
242 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
243
244 //get characteristic length
245 const Scalar characteristicLength = volVars.characteristicLength() ;
246 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
247
248 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
249
250 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
251 {
252 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
253
254 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
255 {
256 if (compIdx <= numPhases)
257 {
258 //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
259 if (phaseIdx == compIdx)
260 continue;
261
262 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
263
264 //additionally get equilibrium values from volume variables
265 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
266 //get the diffusion coefficient
267 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
268
269 //now compute the flux
270 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
271
272 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
273 componentIntoPhaseMassTransfer[compIdx][compIdx] -= compFluxIntoOtherPhase;
274 }
275 }
276 }
277
278 // Actually add the mass transfer to the sources which might
279 // exist externally
280 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
281 {
282 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
283 {
284 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
285 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
286
287 using std::isfinite;
288 if (!isfinite(source[eqIdx]))
289 DUNE_THROW(NumericalProblem, "Calculated non-finite source");
290 }
291 }
292
293 if constexpr (ModelTraits::enableThermalNonEquilibrium())
294 {
295 // Call the (kinetic) Energy module, for the source term.
296 // it has to be called from here, because the mass transferred has to be known.
297 EnergyLocalResidual::computeSourceEnergy(source,
298 element,
299 fvGeometry,
300 elemVolVars,
301 scv);
302 }
303
304 // add contributions from volume flux sources
305 source += problem.source(element, fvGeometry, elemVolVars, scv);
306
307 // add contribution from possible point sources
308 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
309
310 return source;
311 }
312 };
313
314} // end namespace Dumux
315
316#endif
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:62
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:231
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:142
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:177
Definition: porousmediumflow/nonequilibrium/localresidual.hh:26
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:27
Defines all properties used in Dumux.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...