3.4
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>
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, bool enableChemicalNonEquilibrium>
39
40template <class TypeTag>
42
48template<class TypeTag>
49class NonEquilibriumLocalResidualImplementation<TypeTag, false>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
50{
54 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
55 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
57 using Element = typename GridView::template Codim<0>::Entity;
58 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61
62public:
63 using ParentType::ParentType;
64
74 NumEqVector computeSource(const Problem& problem,
75 const Element& element,
76 const FVElementGeometry& fvGeometry,
77 const ElementVolumeVariables& elemVolVars,
78 const SubControlVolume &scv) const
79 {
80 NumEqVector source(0.0);
81
82 // add contributions from volume flux sources
83 source += problem.source(element, fvGeometry, elemVolVars, scv);
84
85 // add contribution from possible point sources
86 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
87
88 // Call the (kinetic) Energy module, for the source term.
89 // it has to be called from here, because the mass transfered has to be known.
90 if constexpr(ModelTraits::enableThermalNonEquilibrium())
91 {
92 EnergyLocalResidual::computeSourceEnergy(source,
93 element,
94 fvGeometry,
95 elemVolVars,
96 scv);
97 }
98
99 return source;
100 }
101
102};
103
108template<class TypeTag>
109class NonEquilibriumLocalResidualImplementation<TypeTag, true>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
110{
114 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
115 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
116 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
119 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
121 using Element = typename GridView::template Codim<0>::Entity;
122 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
126
128 using Indices = typename ModelTraits::Indices;
129
130 static constexpr int numPhases = ModelTraits::numFluidPhases();
131 static constexpr int numComponents = ModelTraits::numFluidComponents();
132
133 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
134 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
135 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
136 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
137
138 static_assert(numPhases > 1,
139 "chemical non-equlibrium only makes sense for multiple phases");
140 static_assert(numPhases == numComponents,
141 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
142 static_assert(ModelTraits::useMoles(),
143 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
144
145public:
146 using ParentType::ParentType;
154 NumEqVector computeStorage(const Problem& problem,
155 const SubControlVolume& scv,
156 const VolumeVariables& volVars) const
157 {
158 NumEqVector storage(0.0);
159
160 // compute storage term of all components within all phases
161 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
162 {
163 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
164 {
165 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
166 storage[eqIdx] += volVars.porosity()
167 * volVars.saturation(phaseIdx)
168 * volVars.molarDensity(phaseIdx)
169 * volVars.moleFraction(phaseIdx, compIdx);
170 }
172 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
173 }
175 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
176 return storage;
177 }
178
189 NumEqVector computeFlux(const Problem& problem,
190 const Element& element,
191 const FVElementGeometry& fvGeometry,
192 const ElementVolumeVariables& elemVolVars,
193 const SubControlVolumeFace& scvf,
194 const ElementFluxVariablesCache& elemFluxVarsCache) const
195 {
196 FluxVariables fluxVars;
197 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
198 NumEqVector flux(0.0);
199
200 // advective fluxes
201 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
202 {
203 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
204 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
205 {
206 // get equation index
207 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
208
209 // the physical quantities for which we perform upwinding
210 const auto upwindTerm = [phaseIdx, compIdx] (const auto& volVars)
211 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
212
213 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
214
215 // do not add diffusive flux of main component, as that is not done in master as well
216 if (compIdx == phaseIdx)
217 continue;
218 //check for the reference system and adapt units of the diffusive flux accordingly.
219 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
220 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
221 else
222 flux[eqIdx] += diffusiveFluxes[compIdx];
223 }
225 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
226 }
227
229 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
230
231 return flux;
232 }
233
243 NumEqVector computeSource(const Problem& problem,
244 const Element& element,
245 const FVElementGeometry& fvGeometry,
246 const ElementVolumeVariables& elemVolVars,
247 const SubControlVolume &scv) const
248 {
249 NumEqVector source(0.0);
250 // In the case of a kinetic consideration, mass transfer
251 // between phases is realized via source terms there is a
252 // balance equation for each component in each phase
253 const auto& volVars = elemVolVars[scv];
254 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
255
256 //get characteristic length
257 const Scalar characteristicLength = volVars.characteristicLength() ;
258 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
259
260 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
261
262 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
263 {
264 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
265
266 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
267 {
268 if (compIdx <= numPhases)
269 {
270 //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
271 if (phaseIdx == compIdx)
272 continue;
273
274 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
275
276 //additionally get equilibrium values from volume variables
277 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
278 //get the diffusion coefficient
279 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, FluidSystem::getMainComponent(phaseIdx), compIdx);
280
281 //now compute the flux
282 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
283
284 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
285 componentIntoPhaseMassTransfer[compIdx][compIdx] -= compFluxIntoOtherPhase;
286 }
287 }
288 }
289
290 // Actually add the mass transfer to the sources which might
291 // exist externally
292 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
293 {
294 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
295 {
296 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
297 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
298
299 using std::isfinite;
300 if (!isfinite(source[eqIdx]))
301 DUNE_THROW(NumericalProblem, "Calculated non-finite source");
302 }
303 }
304
305 if constexpr (ModelTraits::enableThermalNonEquilibrium())
306 {
307 // Call the (kinetic) Energy module, for the source term.
308 // it has to be called from here, because the mass transfered has to be known.
309 EnergyLocalResidual::computeSourceEnergy(source,
310 element,
311 fvGeometry,
312 elemVolVars,
313 scv);
314 }
315
316 // add contributions from volume flux sources
317 source += problem.source(element, fvGeometry, elemVolVars, scv);
318
319 // add contribution from possible point sources
320 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
321
322 return source;
323 }
324 };
325
326} // end namespace Dumux
327
328#endif
A helper to deduce a vector with the same size as numbers of equations.
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:46
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Definition: porousmediumflow/nonequilibrium/localresidual.hh:38
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:74
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:243
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:154
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:189
Declares all properties used in Dumux.
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...