version 3.8
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-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_TRACER_LOCAL_RESIDUAL_HH
15#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
16
17#include <dune/common/exceptions.hh>
18
25
26namespace Dumux {
27
33template<class TypeTag>
34class TracerLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
35{
40 using FVElementGeometry = typename GridGeometry::LocalView;
41 using SubControlVolume = typename GridGeometry::SubControlVolume;
42 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
43 using Extrusion = Extrusion_t<GridGeometry>;
46 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
48 using Element = typename GridView::template Codim<0>::Entity;
49 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
53
54 static constexpr int numComponents = ModelTraits::numFluidComponents();
55 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
56 static constexpr int phaseIdx = 0;
57
58public:
59 using ParentType::ParentType;
60
72 NumEqVector computeStorage(const Problem& problem,
73 const SubControlVolume& scv,
74 const VolumeVariables& volVars) const
75 {
76 NumEqVector storage(0.0);
77
78 // regularize saturation so we don't get singular matrices when the saturation is zero
79 // note that the fluxes will still be zero (zero effective diffusion coefficient),
80 // and we still solve the equation storage = 0 yielding the correct result
81 using std::max;
82 const Scalar saturation = max(1e-8, volVars.saturation(phaseIdx));
83
84 // formulation with mole balances
85 if (useMoles)
86 {
87 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
88 storage[compIdx] += volVars.porosity()
89 * volVars.molarDensity(phaseIdx)
90 * volVars.moleFraction(phaseIdx, compIdx)
91 * saturation;
92 }
93 // formulation with mass balances
94 else
95 {
96 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
97 storage[compIdx] += volVars.porosity()
98 * volVars.density(phaseIdx)
99 * volVars.massFraction(phaseIdx, compIdx)
100 * saturation;
101 }
102
103 return storage;
104 }
105
117 NumEqVector computeFlux(const Problem& problem,
118 const Element& element,
119 const FVElementGeometry& fvGeometry,
120 const ElementVolumeVariables& elemVolVars,
121 const SubControlVolumeFace& scvf,
122 const ElementFluxVariablesCache& elemFluxVarsCache) const
123 {
124 FluxVariables fluxVars;
125 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
126
127 NumEqVector flux(0.0);
128 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
129
130 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
131
132 if (useMoles) // formulation with mole balances
133 {
134 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
135 {
136 // the physical quantities for which we perform upwinding
137 auto upwindTerm = [compIdx](const auto& volVars)
138 { return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
139
140 // advective fluxes
141 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
142 // diffusive fluxes
143 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
144 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
145 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
146 flux[compIdx] += diffusiveFluxes[compIdx];
147 else
148 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
149 }
150 }
151 else // formulation with mass balances
152 {
153 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
154 {
155 // the physical quantities for which we perform upwinding
156 auto upwindTerm = [compIdx](const auto& volVars)
157 { return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
158
159 // advective fluxes
160 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
161 // diffusive fluxes
162 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
163 flux[compIdx] += diffusiveFluxes[compIdx];
164 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
165 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
166 else
167 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
168 }
169 }
170
171 if constexpr (ModelTraits::enableCompositionalDispersion())
172 {
173 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
174 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
175 {
176 flux[compIdx] += dispersionFluxes[compIdx];
177 }
178 }
179
180 return flux;
181 }
182
193 template<class PartialDerivativeMatrix>
194 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
195 const Problem& problem,
196 const Element& element,
197 const FVElementGeometry& fvGeometry,
198 const VolumeVariables& curVolVars,
199 const SubControlVolume& scv) const
200 {
201 // regularize saturation so we don't get singular matrices when the saturation is zero
202 // note that the fluxes will still be zero (zero effective diffusion coefficient),
203 // and we still solve the equation storage = 0 yielding the correct result
204 using std::max;
205 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
206
207 const auto porosity = curVolVars.porosity();
208 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
209 const auto d_storage = Extrusion::volume(fvGeometry, scv)*porosity*rho*saturation/this->timeLoop().timeStepSize();
210
211 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
212 partialDerivatives[compIdx][compIdx] += d_storage;
213 }
214
225 template<class PartialDerivativeMatrix>
226 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
227 const Problem& problem,
228 const Element& element,
229 const FVElementGeometry& fvGeometry,
230 const VolumeVariables& curVolVars,
231 const SubControlVolume& scv) const
232 {
233 // TODO maybe forward to the problem? -> necessary for reaction terms
234 }
235
236 template<class PartialDerivativeMatrices, class T = TypeTag>
237 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethods::box, void>
238 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
239 const Problem& problem,
240 const Element& element,
241 const FVElementGeometry& fvGeometry,
242 const ElementVolumeVariables& curElemVolVars,
243 const ElementFluxVariablesCache& elemFluxVarsCache,
244 const SubControlVolumeFace& scvf) const
245 {
246 if constexpr (FVElementGeometry::GridGeometry::discMethod != DiscretizationMethods::cctpfa)
247 DUNE_THROW(Dune::NotImplemented, "Analytic flux differentiation only implemented for tpfa");
248
249 // advective term: we do the same for all tracer components
250 auto rho = [](const VolumeVariables& volVars)
251 { return useMoles ? volVars.molarDensity() : volVars.density(); };
252
253 // the volume flux
254 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
255
256 // the upwind weight
257 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
258
259 // get the inside and outside volvars
260 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
261 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
262
263 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
264 const Scalar outsideWeight = 1.0 - insideWeight;
265 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
266 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
267
268 // diffusive term
269 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
270 const auto& fluxCache = elemFluxVarsCache[scvf];
271 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
272 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
273 const Scalar massOrMolarDensity = 0.5*(rhoInside + rhoOutside);
274
275 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
276 {
277 // diffusive term
278 Scalar diffDeriv = 0.0;
279 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
280 {
281 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
282 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx);
283 }
284 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
285 {
286 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)
287 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
288 }
289 else
290 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
291
292 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
293 if (!scvf.boundary())
294 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
295 }
296 }
297
298 template<class JacobianMatrix, class T = TypeTag>
299 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
300 addFluxDerivatives(JacobianMatrix& A,
301 const Problem& problem,
302 const Element& element,
303 const FVElementGeometry& fvGeometry,
304 const ElementVolumeVariables& curElemVolVars,
305 const ElementFluxVariablesCache& elemFluxVarsCache,
306 const SubControlVolumeFace& scvf) const
307 {
308
309 // advective term: we do the same for all tracer components
310 auto rho = [](const VolumeVariables& volVars)
311 { return useMoles ? volVars.molarDensity() : volVars.density(); };
312
313 // the volume flux
314 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
315
316 // the upwind weight
317 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
318
319 // get the inside and outside volvars
320 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
321 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
322
323 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
324 const auto outsideWeight = 1.0 - insideWeight;
325 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
326 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
327
328 // diffusive term
329 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
331 const auto ti = DiffusionType::calculateTransmissibilities(problem,
332 element,
333 fvGeometry,
334 curElemVolVars,
335 scvf,
336 elemFluxVarsCache[scvf],
337 phaseIdx);
338 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
339 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
340
341 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
342 {
343 for (const auto& scv : scvs(fvGeometry))
344 {
345 // diffusive term
346 auto diffDeriv = 0.0;
347 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
348 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
349 : ti[compIdx][scv.indexInElement()];
350 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
351 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
352 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
353 else
354 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
355 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
356 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
357 }
358
359 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
360 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
361 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
362 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
363 }
364 }
365
366 template<class PartialDerivativeMatrices>
367 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
368 const Problem& problem,
369 const Element& element,
370 const FVElementGeometry& fvGeometry,
371 const ElementVolumeVariables& curElemVolVars,
372 const ElementFluxVariablesCache& elemFluxVarsCache,
373 const SubControlVolumeFace& scvf) const
374 {
375 // do the same as for inner facets
376 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry,
377 curElemVolVars, elemFluxVarsCache, scvf);
378 }
379
380 template<class PartialDerivativeMatrices>
381 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
382 const Problem& problem,
383 const Element& element,
384 const FVElementGeometry& fvGeometry,
385 const ElementVolumeVariables& curElemVolVars,
386 const ElementFluxVariablesCache& elemFluxVarsCache,
387 const SubControlVolumeFace& scvf) const
388 {
389 // TODO maybe forward to the problem?
390 }
391};
392
393} // end namespace Dumux
394
395#endif
Element-wise calculation of the local residual for problems using fully implicit tracer model.
Definition: porousmediumflow/tracer/localresidual.hh:35
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:72
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:381
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:238
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:226
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:194
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:117
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:300
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:367
Defines all properties used in Dumux.
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 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
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
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
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:31
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:127
Definition: adapt.hh:17
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...