version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
porousmediumflow/compositional/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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
15#define DUMUX_COMPOSITIONAL_LOCAL_RESIDUAL_HH
16
17#include <vector>
18#include <dune/common/exceptions.hh>
24
25namespace Dumux {
26
32template<class TypeTag>
35{
40 using FVElementGeometry = typename GridGeometry::LocalView;
41 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
42 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
45 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
46 using GridView = typename GridGeometry::GridView;
47 using Element = typename GridView::template Codim<0>::Entity;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
53 using Indices = typename ModelTraits::Indices;
54
55 static constexpr int numPhases = ModelTraits::numFluidPhases();
56 static constexpr int numComponents = ModelTraits::numFluidComponents();
57 static constexpr bool useMoles = ModelTraits::useMoles();
58
59 enum { conti0EqIdx = Indices::conti0EqIdx };
60
62 static constexpr int replaceCompEqIdx = ModelTraits::replaceCompEqIdx();
63 static constexpr bool useTotalMoleOrMassBalance = replaceCompEqIdx < numComponents;
64
65public:
66 using ParentType::ParentType;
67
79 NumEqVector computeStorage(const Problem& problem,
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars) const
82 {
83 NumEqVector storage(0.0);
84
85 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
86 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
87
88 const auto massOrMoleFraction= [](const auto& volVars, const int phaseIdx, const int compIdx)
89 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
90
91 // compute storage term of all components within all phases
92 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
93 {
94 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
95 {
96 auto eqIdx = conti0EqIdx + compIdx;
97 if (eqIdx != replaceCompEqIdx)
98 storage[eqIdx] += volVars.porosity()
99 * volVars.saturation(phaseIdx)
100 * massOrMoleDensity(volVars, phaseIdx)
101 * massOrMoleFraction(volVars, phaseIdx, compIdx);
102 }
103
104 // in case one balance is substituted by the total mole balance
105 if (useTotalMoleOrMassBalance)
106 storage[replaceCompEqIdx] += massOrMoleDensity(volVars, phaseIdx)
107 * volVars.porosity()
108 * volVars.saturation(phaseIdx);
109
111 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, phaseIdx);
112 }
113
115 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
116
117 return storage;
118 }
119
131 NumEqVector computeFlux(const Problem& problem,
132 const Element& element,
133 const FVElementGeometry& fvGeometry,
134 const ElementVolumeVariables& elemVolVars,
135 const SubControlVolumeFace& scvf,
136 const ElementFluxVariablesCache& elemFluxVarsCache) const
137 {
138 FluxVariables fluxVars;
139 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
140 static constexpr auto referenceSystemFormulationDiffusion = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
141 // get upwind weights into local scope
142 NumEqVector flux(0.0);
143
144 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
145 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
146
147 const auto massOrMoleFraction = [](const auto& volVars, const int phaseIdx, const int compIdx)
148 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
149
150 // Check for the reference system and adapt units of the flux accordingly.
151 const auto adaptFluxUnits = [](const Scalar referenceFlux, const Scalar molarMass,
152 const ReferenceSystemFormulation referenceSystemFormulation)
153 {
154 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
155 return useMoles ? referenceFlux/molarMass
156 : referenceFlux;
157 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
158 return useMoles ? referenceFlux
159 : referenceFlux*molarMass;
160 else
161 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
162 };
163
164 for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
165 {
166 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
167 auto dispersiveFluxes = decltype(diffusiveFluxes)(0.0);
168 if constexpr (ModelTraits::enableCompositionalDispersion())
169 dispersiveFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
170
171 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
172 {
173 // get equation index
174 const auto eqIdx = conti0EqIdx + compIdx;
175
176 // the physical quantities for which we perform upwinding
177 const auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, phaseIdx, compIdx] (const auto& volVars)
178 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx)*volVars.mobility(phaseIdx); };
179
180 // advective fluxes (only for the component balances)
181 if (eqIdx != replaceCompEqIdx)
182 flux[eqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
183
184 // diffusive and dispersive fluxes
185 auto diffusiveAndDispersiveFluxes = adaptFluxUnits(diffusiveFluxes[compIdx],
186 FluidSystem::molarMass(compIdx),
187 referenceSystemFormulationDiffusion);
188 if constexpr (ModelTraits::enableCompositionalDispersion())
189 {
190 static constexpr auto referenceSystemFormulationDispersion =
191 FluxVariables::DispersionFluxType::referenceSystemFormulation();
192 diffusiveAndDispersiveFluxes += adaptFluxUnits(dispersiveFluxes[compIdx],
193 FluidSystem::molarMass(compIdx),
194 referenceSystemFormulationDispersion);
195 }
196 if(eqIdx != replaceCompEqIdx)
197 flux[eqIdx] += diffusiveAndDispersiveFluxes;
198 if (useTotalMoleOrMassBalance)
199 flux[replaceCompEqIdx] += diffusiveAndDispersiveFluxes;
200 }
201
202 // in case one balance is substituted by the total mole balance
203 if (useTotalMoleOrMassBalance)
204 {
205 // the physical quantities for which we perform upwinding
206 const auto upwindTerm = [&massOrMoleDensity, phaseIdx] (const auto& volVars)
207 { return massOrMoleDensity(volVars, phaseIdx)*volVars.mobility(phaseIdx); };
208
209 flux[replaceCompEqIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
210 }
211
213 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, phaseIdx);
214 }
215
217 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
218 EnergyLocalResidual::heatDispersionFlux(flux, fluxVars);
219
220 return flux;
221 }
222};
223
224} // end namespace Dumux
225
226#endif
Element-wise calculation of the local residual for problems using compositional fully implicit model.
Definition: porousmediumflow/compositional/localresidual.hh:35
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluates the total flux of all conservation quantities over a face of a sub-control volume.
Definition: porousmediumflow/compositional/localresidual.hh:131
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluates the amount of all conservation quantities (e.g. phase mass) within a sub-control volume.
Definition: porousmediumflow/compositional/localresidual.hh:79
Defines all properties used in Dumux.
The default local operator than can be specialized for each discretization scheme.
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
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:54
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition: defaultlocaloperator.hh:27
A helper to deduce a vector with the same size as numbers of equations.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...