version 3.10.0
Loading...
Searching...
No Matches
porousmediumflow/tracer/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//
13
14#ifndef DUMUX_TRACER_LOCAL_RESIDUAL_HH
15#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
16
17#include <dune/common/exceptions.hh>
18
26
27namespace Dumux {
28
34template<class TypeTag>
37{
42 using FVElementGeometry = typename GridGeometry::LocalView;
43 using SubControlVolume = typename GridGeometry::SubControlVolume;
44 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
45 using Extrusion = Extrusion_t<GridGeometry>;
48 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
49 using GridView = typename GridGeometry::GridView;
50 using Element = typename GridView::template Codim<0>::Entity;
51 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
55
56 static constexpr int numComponents = ModelTraits::numFluidComponents();
57 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
58 static constexpr int phaseIdx = 0;
59
60public:
61 using ParentType::ParentType;
62
74 NumEqVector computeStorage(const Problem& problem,
75 const SubControlVolume& scv,
76 const VolumeVariables& volVars) const
77 {
78 NumEqVector storage(0.0);
79
80 // regularize saturation so we don't get singular matrices when the saturation is zero
81 // note that the fluxes will still be zero (zero effective diffusion coefficient),
82 // and we still solve the equation storage = 0 yielding the correct result
83 using std::max;
84 const Scalar saturation = max(1e-8, volVars.saturation(phaseIdx));
85
86 // formulation with mole balances
87 if (useMoles)
88 {
89 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
90 storage[compIdx] += volVars.porosity()
91 * volVars.molarDensity(phaseIdx)
92 * volVars.moleFraction(phaseIdx, compIdx)
93 * saturation;
94 }
95 // formulation with mass balances
96 else
97 {
98 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
99 storage[compIdx] += volVars.porosity()
100 * volVars.density(phaseIdx)
101 * volVars.massFraction(phaseIdx, compIdx)
102 * saturation;
103 }
104
105 return storage;
106 }
107
119 NumEqVector computeFlux(const Problem& problem,
120 const Element& element,
121 const FVElementGeometry& fvGeometry,
122 const ElementVolumeVariables& elemVolVars,
123 const SubControlVolumeFace& scvf,
124 const ElementFluxVariablesCache& elemFluxVarsCache) const
125 {
126 FluxVariables fluxVars;
127 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
128 static constexpr auto referenceSystemFormulationDiffusion = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
129
130 NumEqVector flux(0.0);
131
132 const auto massOrMoleDensity = [](const auto& volVars, const int phaseIdx)
133 { return useMoles ? volVars.molarDensity(phaseIdx) : volVars.density(phaseIdx); };
134
135 const auto massOrMoleFraction = [](const auto& volVars, const int phaseIdx, const int compIdx)
136 { return useMoles ? volVars.moleFraction(phaseIdx, compIdx) : volVars.massFraction(phaseIdx, compIdx); };
137
138 // Check for the reference system and adapt units of the flux accordingly.
139 const auto adaptFluxUnits = [](const Scalar referenceFlux, const Scalar molarMass,
140 const ReferenceSystemFormulation referenceSystemFormulation)
141 {
142 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
143 return useMoles ? referenceFlux/molarMass
144 : referenceFlux;
145 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
146 return useMoles ? referenceFlux
147 : referenceFlux*molarMass;
148 else
149 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
150 };
151
152 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
153 auto dispersiveFluxes = decltype(diffusiveFluxes)(0.0);
154 if constexpr (ModelTraits::enableCompositionalDispersion())
155 dispersiveFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
156
157 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
158 {
159 // the physical quantities for which we perform upwinding
160 auto upwindTerm = [&massOrMoleDensity, &massOrMoleFraction, compIdx](const auto& volVars)
161 { return massOrMoleDensity(volVars, phaseIdx)*massOrMoleFraction(volVars, phaseIdx, compIdx); };
162
163 // advective fluxes
164 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
165 // diffusive and dispersive fluxes
166 flux[compIdx] += adaptFluxUnits(diffusiveFluxes[compIdx],
167 FluidSystem::molarMass(compIdx),
168 referenceSystemFormulationDiffusion);
169 if constexpr (ModelTraits::enableCompositionalDispersion())
170 {
171 static constexpr auto referenceSystemFormulationDispersion =
172 FluxVariables::DispersionFluxType::referenceSystemFormulation();
173 flux[compIdx] += adaptFluxUnits(dispersiveFluxes[compIdx],
174 FluidSystem::molarMass(compIdx),
175 referenceSystemFormulationDispersion);
176 }
177 }
178
179 return flux;
180 }
181
192 template<class PartialDerivativeMatrix>
193 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
194 const Problem& problem,
195 const Element& element,
196 const FVElementGeometry& fvGeometry,
197 const VolumeVariables& curVolVars,
198 const SubControlVolume& scv) const
199 {
200 // regularize saturation so we don't get singular matrices when the saturation is zero
201 // note that the fluxes will still be zero (zero effective diffusion coefficient),
202 // and we still solve the equation storage = 0 yielding the correct result
203 using std::max;
204 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
205
206 const auto porosity = curVolVars.porosity();
207 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
208 const auto d_storage = Extrusion::volume(fvGeometry, scv)*porosity*rho*saturation/this->timeLoop().timeStepSize();
209
210 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
211 partialDerivatives[compIdx][compIdx] += d_storage;
212 }
213
224 template<class PartialDerivativeMatrix>
225 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
226 const Problem& problem,
227 const Element& element,
228 const FVElementGeometry& fvGeometry,
229 const VolumeVariables& curVolVars,
230 const SubControlVolume& scv) const
231 {
232 // TODO maybe forward to the problem? -> necessary for reaction terms
233 }
234
235 template<class PartialDerivativeMatrices, class T = TypeTag>
236 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethods::box, void>
237 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
238 const Problem& problem,
239 const Element& element,
240 const FVElementGeometry& fvGeometry,
241 const ElementVolumeVariables& curElemVolVars,
242 const ElementFluxVariablesCache& elemFluxVarsCache,
243 const SubControlVolumeFace& scvf) const
244 {
245 if constexpr (FVElementGeometry::GridGeometry::discMethod != DiscretizationMethods::cctpfa)
246 DUNE_THROW(Dune::NotImplemented, "Analytic flux differentiation only implemented for tpfa");
247
248 // advective term: we do the same for all tracer components
249 auto rho = [](const VolumeVariables& volVars)
250 { return useMoles ? volVars.molarDensity() : volVars.density(); };
251
252 // the volume flux
253 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
254
255 // the upwind weight
256 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
257
258 // get the inside and outside volvars
259 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
260 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
261
262 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
263 const Scalar outsideWeight = 1.0 - insideWeight;
264 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
265 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
266
267 // diffusive term
268 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
269 const auto& fluxCache = elemFluxVarsCache[scvf];
270 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
271 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
272 const Scalar massOrMolarDensity = 0.5*(rhoInside + rhoOutside);
273
274 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
275 {
276 // diffusive term
277 Scalar diffDeriv = 0.0;
278 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
279 {
280 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
281 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx);
282 }
283 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
284 {
285 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)
286 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
287 }
288 else
289 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
290
291 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
292 if (!scvf.boundary())
293 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
294 }
295 }
296
297 template<class JacobianMatrix, class T = TypeTag>
298 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
299 addFluxDerivatives(JacobianMatrix& A,
300 const Problem& problem,
301 const Element& element,
302 const FVElementGeometry& fvGeometry,
303 const ElementVolumeVariables& curElemVolVars,
304 const ElementFluxVariablesCache& elemFluxVarsCache,
305 const SubControlVolumeFace& scvf) const
306 {
307
308 // advective term: we do the same for all tracer components
309 auto rho = [](const VolumeVariables& volVars)
310 { return useMoles ? volVars.molarDensity() : volVars.density(); };
311
312 // the volume flux
313 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
314
315 // the upwind weight
316 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
317
318 // get the inside and outside volvars
319 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
320 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
321
322 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
323 const auto outsideWeight = 1.0 - insideWeight;
324 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
325 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
326
327 // diffusive term
328 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
330 const auto ti = DiffusionType::calculateTransmissibilities(problem,
331 element,
332 fvGeometry,
333 curElemVolVars,
334 scvf,
335 elemFluxVarsCache[scvf],
336 phaseIdx);
337 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
338 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
339
340 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
341 {
342 for (const auto& scv : scvs(fvGeometry))
343 {
344 // diffusive term
345 auto diffDeriv = 0.0;
346 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
347 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
348 : ti[compIdx][scv.indexInElement()];
349 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
350 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
351 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
352 else
353 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
354 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
355 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
356 }
357
358 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
359 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
360 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
361 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
362 }
363 }
364
365 template<class PartialDerivativeMatrices>
366 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
367 const Problem& problem,
368 const Element& element,
369 const FVElementGeometry& fvGeometry,
370 const ElementVolumeVariables& curElemVolVars,
371 const ElementFluxVariablesCache& elemFluxVarsCache,
372 const SubControlVolumeFace& scvf) const
373 {
374 // do the same as for inner facets
375 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry,
376 curElemVolVars, elemFluxVarsCache, scvf);
377 }
378
379 template<class PartialDerivativeMatrices>
380 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
381 const Problem& problem,
382 const Element& element,
383 const FVElementGeometry& fvGeometry,
384 const ElementVolumeVariables& curElemVolVars,
385 const ElementFluxVariablesCache& elemFluxVarsCache,
386 const SubControlVolumeFace& scvf) const
387 {
388 // TODO maybe forward to the problem?
389 }
390};
391
392} // end namespace Dumux
393
394#endif
Element-wise calculation of the local residual for problems using fully implicit tracer model.
Definition porousmediumflow/tracer/localresidual.hh:37
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/tracer/localresidual.hh:74
void addRobinFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition porousmediumflow/tracer/localresidual.hh:380
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod !=DiscretizationMethods::box, void > addFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition porousmediumflow/tracer/localresidual.hh:237
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
TODO docme!
Definition porousmediumflow/tracer/localresidual.hh:225
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
TODO docme!
Definition porousmediumflow/tracer/localresidual.hh:193
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/tracer/localresidual.hh:119
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::box, void > addFluxDerivatives(JacobianMatrix &A, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition porousmediumflow/tracer/localresidual.hh:299
void addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition porousmediumflow/tracer/localresidual.hh:366
Defines all properties used in Dumux.
The default local operator than can be specialized for each discretization scheme.
Helper classes to compute the integration elements.
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
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition referencesystemformulation.hh:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition referencesystemformulation.hh:33
@ massAveraged
Definition referencesystemformulation.hh:34
@ molarAveraged
Definition referencesystemformulation.hh:34
T getParamFromGroup(Args &&... args)
A free function to get a parameter from the parameter tree singleton with a model group.
Definition parameters.hh:149
T getParam(Args &&... args)
A free function to get a parameter from the parameter tree singleton.
Definition parameters.hh:139
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:310
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition propertysystem.hh:296
The available discretization methods in Dumux.
constexpr CCTpfa cctpfa
Definition method.hh:145
constexpr Box box
Definition method.hh:147
Definition adapt.hh:17
typename Detail::DiscretizationDefaultLocalOperator< TypeTag >::type DiscretizationDefaultLocalOperator
Definition defaultlocaloperator.hh:26
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition extrusion.hh:166
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...