3.1-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
31
32namespace Dumux {
33
34// forward declaration
35template<class TypeTag, bool enableChemicalNonEquilibrium>
37
38template <class TypeTag>
40
45template<class TypeTag>
46class NonEquilibriumLocalResidualImplementation<TypeTag, false>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
47{
52 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
53 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
54 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
56 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
58 using Element = typename GridView::template Codim<0>::Entity;
59 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
63 using Indices = typename ModelTraits::Indices;
64
65 static constexpr int numPhases = ModelTraits::numFluidPhases();
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67 enum { conti0EqIdx = Indices::conti0EqIdx };
68public:
69 using ParentType::ParentType;
70
86 NumEqVector computeFlux(const Problem& problem,
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf,
91 const ElementFluxVariablesCache& elemFluxVarsCache) const
92 {
93 FluxVariables fluxVars;
94 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
95 // get upwind weights into local scope
96 NumEqVector flux(0.0);
97
98 const auto moleDensity = [](const auto& volVars, const int phaseIdx)
99 { return volVars.molarDensity(phaseIdx); };
100
101 const auto moleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
102 { return volVars.moleFraction(phaseIdx, compIdx); };
103
104 // advective fluxes
105 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
106 {
107 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
108 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
109 {
110 // get equation index
111 const auto eqIdx = conti0EqIdx + compIdx;
112
113 // the physical quantities for which we perform upwinding
114 const auto upwindTerm = [&moleDensity, &moleFraction, phaseIdx, compIdx] (const auto& volVars)
115 { return moleDensity(volVars, phaseIdx)*moleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
116
117 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
118
119 // diffusive fluxes (only for the component balances)
120 //check for the reference system and adapt units of the diffusive flux accordingly.
121 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
122 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
123 else
124 flux[eqIdx] += diffusiveFluxes[compIdx];
125 }
126
128 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
129 }
131 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
132
133 return flux;
134
135
136 }
137
138 NumEqVector computeSource(const Problem& problem,
139 const Element& element,
140 const FVElementGeometry& fvGeometry,
141 const ElementVolumeVariables& elemVolVars,
142 const SubControlVolume &scv) const
143 {
144 NumEqVector source(0.0);
145
146 // add contributions from volume flux sources
147 source += problem.source(element, fvGeometry, elemVolVars, scv);
148
149 // add contribution from possible point sources
150 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
151 // Call the (kinetic) Energy module, for the source term.
152 // it has to be called from here, because the mass transfered has to be known.
153 EnergyLocalResidual::computeSourceEnergy(source,
154 element,
155 fvGeometry,
156 elemVolVars,
157 scv);
158 return source;
159 }
160
161};
162
166template<class TypeTag>
167class NonEquilibriumLocalResidualImplementation<TypeTag, true>: public GetPropType<TypeTag, Properties::EquilibriumLocalResidual>
168{
173 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
174 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
175 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
178 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
180 using Element = typename GridView::template Codim<0>::Entity;
181 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
186
188 using Indices = typename ModelTraits::Indices;
189
190 static constexpr int numPhases = ModelTraits::numFluidPhases();
191 static constexpr int numComponents = ModelTraits::numFluidComponents();
192 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
193
194 static constexpr auto conti0EqIdx = Indices::conti0EqIdx;
195 static constexpr auto comp1Idx = FluidSystem::comp1Idx;
196 static constexpr auto comp0Idx = FluidSystem::comp0Idx;
197 static constexpr auto phase0Idx = FluidSystem::phase0Idx;
198 static constexpr auto phase1Idx = FluidSystem::phase1Idx;
199
200 static_assert(numPhases == numComponents,
201 "currently chemical non-equilibrium is only available when numPhases equals numComponents");
202 static_assert(useMoles == true,
203 "chemical nonequilibrium can only be calculated based on mole fractions not mass fractions");
204
205public:
206 using ParentType::ParentType;
214 NumEqVector computeStorage(const Problem& problem,
215 const SubControlVolume& scv,
216 const VolumeVariables& volVars) const
217 {
218 NumEqVector storage(0.0);
219
220 // compute storage term of all components within all phases
221 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
222 {
223 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
224 {
225 auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
226 storage[eqIdx] += volVars.porosity()
227 * volVars.saturation(phaseIdx)
228 * volVars.molarDensity(phaseIdx)
229 * volVars.moleFraction(phaseIdx, compIdx);
230 }
232 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, phaseIdx);
233 }
235 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
236 return storage;
237 }
238
249 NumEqVector computeFlux(const Problem& problem,
250 const Element& element,
251 const FVElementGeometry& fvGeometry,
252 const ElementVolumeVariables& elemVolVars,
253 const SubControlVolumeFace& scvf,
254 const ElementFluxVariablesCache& elemFluxVarsCache) const
255 {
256 FluxVariables fluxVars;
257 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
258 // get upwind weights into local scope
259 NumEqVector flux(0.0);
260
261 // advective fluxes
262 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
263 {
264 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
265 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
266 {
267 // get equation index
268 const auto eqIdx = conti0EqIdx + phaseIdx*numComponents + compIdx;
269
270 // the physical quantities for which we perform upwinding
271 const auto upwindTerm = [phaseIdx, compIdx] (const auto& volVars)
272 { return volVars.molarDensity(phaseIdx)*volVars.moleFraction(phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
273
274 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
275
276 // do not add diffusive flux of main component, as that is not done in master as well
277 if (compIdx == phaseIdx)
278 continue;
279 //check for the reference system and adapt units of the diffusive flux accordingly.
280 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
281 flux[eqIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
282 else
283 flux[eqIdx] += diffusiveFluxes[compIdx];
284 }
286 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
287 }
288
290 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
291
292 return flux;
293 }
294
304 NumEqVector computeSource(const Problem& problem,
305 const Element& element,
306 const FVElementGeometry& fvGeometry,
307 const ElementVolumeVariables& elemVolVars,
308 const SubControlVolume &scv) const
309 {
310 NumEqVector source(0.0);
311 // In the case of a kinetic consideration, mass transfer
312 // between phases is realized via source terms there is a
313 // balance equation for each component in each phase
314 const auto& volVars = elemVolVars[scv];
315 std::array<std::array<Scalar, numComponents>, numPhases> componentIntoPhaseMassTransfer = {{{0.0},{0.0}}};
316
317 //get characteristic length
318 const Scalar characteristicLength = volVars.characteristicLength() ;
319 const Scalar factorMassTransfer = volVars.factorMassTransfer() ;
320
321 const Scalar awn = volVars.interfacialArea(phase0Idx, phase1Idx);
322
323 for(int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
324 {
325 const Scalar sherwoodNumber = volVars.sherwoodNumber(phaseIdx);
326
327 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
328 {
329 if (compIdx <= numPhases)
330 {
331 //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
332 if (phaseIdx == compIdx)
333 continue;
334
335 const Scalar xNonEquil = volVars.moleFraction(phaseIdx, compIdx);
336
337 //additionally get equilibrium values from volume variables
338 const Scalar xEquil = volVars.xEquil(phaseIdx, compIdx);
339 //get the diffusion coefficient
340 const Scalar diffCoeff = volVars.diffusionCoefficient(phaseIdx, compIdx);
341
342 //now compute the flux
343 const Scalar compFluxIntoOtherPhase = factorMassTransfer * (xEquil-xNonEquil)/characteristicLength * awn * volVars.molarDensity(phaseIdx) * diffCoeff * sherwoodNumber;
344
345 componentIntoPhaseMassTransfer[phaseIdx][compIdx] += compFluxIntoOtherPhase;
346 componentIntoPhaseMassTransfer[compIdx][compIdx] += -compFluxIntoOtherPhase;
347 }
348 }
349 }
350
351 // Actually add the mass transfer to the sources which might
352 // exist externally
353 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
354 {
355 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
356 {
357 const unsigned int eqIdx = conti0EqIdx + compIdx + phaseIdx*numComponents;
358 source[eqIdx] += componentIntoPhaseMassTransfer[phaseIdx][compIdx];
359
360 using std::isfinite;
361 if (!isfinite(source[eqIdx]))
362 DUNE_THROW(NumericalProblem, "Calculated non-finite source");
363 }
364 }
365
366 if (ModelTraits::enableThermalNonEquilibrium())
367 {
368 // Call the (kinetic) Energy module, for the source term.
369 // it has to be called from here, because the mass transfered has to be known.
370 EnergyLocalResidual::computeSourceEnergy(source,
371 element,
372 fvGeometry,
373 elemVolVars,
374 scv);
375 }
376
377 // add contributions from volume flux sources
378 source += problem.source(element, fvGeometry, elemVolVars, scv);
379
380 // add contribution from possible point sources
381 source += problem.scvPointSources(element, fvGeometry, elemVolVars, scv);
383
384 return source;
385 }
386 };
387
388} // end namespace Dumux
389
390#endif
bool CheckDefined(const T &value)
Make valgrind complain if the object occupied by an object is undefined.
Definition: valgrind.hh:72
make the local view function available whenever we use the grid geometry
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
std::string moleFraction(int phaseIdx, int compIdx) noexcept
I/O name of mole fraction.
Definition: name.hh:110
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Definition: porousmediumflow/nonequilibrium/localresidual.hh:36
NumEqVector computeSource(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Definition: porousmediumflow/nonequilibrium/localresidual.hh:138
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Calculates the source term of the equation.
Definition: porousmediumflow/nonequilibrium/localresidual.hh:86
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:304
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:214
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:249
Declares all properties used in Dumux.
This file contains the parts of the local residual to calculate the heat conservation in the thermal ...