3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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_TRACER_LOCAL_RESIDUAL_HH
27#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
28
29#include <dune/common/exceptions.hh>
30
38
39namespace Dumux {
40
46template<class TypeTag>
47class TracerLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
48{
53 using FVElementGeometry = typename GridGeometry::LocalView;
54 using SubControlVolume = typename GridGeometry::SubControlVolume;
55 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
56 using Extrusion = Extrusion_t<GridGeometry>;
59 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
61 using Element = typename GridView::template Codim<0>::Entity;
62 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
66
67 static constexpr int numComponents = ModelTraits::numFluidComponents();
68 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
69 static constexpr int phaseIdx = 0;
70
71public:
72 using ParentType::ParentType;
73
85 NumEqVector computeStorage(const Problem& problem,
86 const SubControlVolume& scv,
87 const VolumeVariables& volVars) const
88 {
89 NumEqVector storage(0.0);
90
91 // regularize saturation so we don't get singular matrices when the saturation is zero
92 // note that the fluxes will still be zero (zero effective diffusion coefficient),
93 // and we still solve the equation storage = 0 yielding the correct result
94 using std::max;
95 const Scalar saturation = max(1e-8, volVars.saturation(phaseIdx));
96
97 // formulation with mole balances
98 if (useMoles)
99 {
100 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
101 storage[compIdx] += volVars.porosity()
102 * volVars.molarDensity(phaseIdx)
103 * volVars.moleFraction(phaseIdx, compIdx)
104 * saturation;
105 }
106 // formulation with mass balances
107 else
108 {
109 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
110 storage[compIdx] += volVars.porosity()
111 * volVars.density(phaseIdx)
112 * volVars.massFraction(phaseIdx, compIdx)
113 * saturation;
114 }
115
116 return storage;
117 }
118
130 NumEqVector computeFlux(const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& elemVolVars,
134 const SubControlVolumeFace& scvf,
135 const ElementFluxVariablesCache& elemFluxVarsCache) const
136 {
137 FluxVariables fluxVars;
138 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
139
140 NumEqVector flux(0.0);
141 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
142
143 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
144
145 if (useMoles) // formulation with mole balances
146 {
147 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
148 {
149 // the physical quantities for which we perform upwinding
150 auto upwindTerm = [compIdx](const auto& volVars)
151 { return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
152
153 // advective fluxes
154 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
155 // diffusive fluxes
156 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
157 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
158 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
159 flux[compIdx] += diffusiveFluxes[compIdx];
160 else
161 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
162 }
163 }
164 else // formulation with mass balances
165 {
166 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
167 {
168 // the physical quantities for which we perform upwinding
169 auto upwindTerm = [compIdx](const auto& volVars)
170 { return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
171
172 // advective fluxes
173 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
174 // diffusive fluxes
175 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
176 flux[compIdx] += diffusiveFluxes[compIdx];
177 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
178 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
179 else
180 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
181 }
182 }
183
184 if constexpr (Deprecated::hasEnableCompositionalDispersion<ModelTraits>())
185 {
186 if constexpr (ModelTraits::enableCompositionalDispersion())
187 {
188 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
189 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
190 {
191 flux[compIdx] += dispersionFluxes[compIdx];
192 }
193 }
194 }
195 else
196 enableCompositionalDispersionMissing_<ModelTraits>();
197
198 return flux;
199 }
200
211 template<class PartialDerivativeMatrix>
212 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
213 const Problem& problem,
214 const Element& element,
215 const FVElementGeometry& fvGeometry,
216 const VolumeVariables& curVolVars,
217 const SubControlVolume& scv) const
218 {
219 // regularize saturation so we don't get singular matrices when the saturation is zero
220 // note that the fluxes will still be zero (zero effective diffusion coefficient),
221 // and we still solve the equation storage = 0 yielding the correct result
222 using std::max;
223 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
224
225 const auto porosity = curVolVars.porosity();
226 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
227 const auto d_storage = Extrusion::volume(scv)*porosity*rho*saturation/this->timeLoop().timeStepSize();
228
229 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
230 partialDerivatives[compIdx][compIdx] += d_storage;
231 }
232
243 template<class PartialDerivativeMatrix>
244 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
245 const Problem& problem,
246 const Element& element,
247 const FVElementGeometry& fvGeometry,
248 const VolumeVariables& curVolVars,
249 const SubControlVolume& scv) const
250 {
251 // TODO maybe forward to the problem? -> necessary for reaction terms
252 }
253
254 template<class PartialDerivativeMatrices, class T = TypeTag>
255 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethods::box, void>
256 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
257 const Problem& problem,
258 const Element& element,
259 const FVElementGeometry& fvGeometry,
260 const ElementVolumeVariables& curElemVolVars,
261 const ElementFluxVariablesCache& elemFluxVarsCache,
262 const SubControlVolumeFace& scvf) const
263 {
264 if constexpr (FVElementGeometry::GridGeometry::discMethod != DiscretizationMethods::cctpfa)
265 DUNE_THROW(Dune::NotImplemented, "Analytic flux differentiation only implemented for tpfa");
266
267 // advective term: we do the same for all tracer components
268 auto rho = [](const VolumeVariables& volVars)
269 { return useMoles ? volVars.molarDensity() : volVars.density(); };
270
271 // the volume flux
272 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
273
274 // the upwind weight
275 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
276
277 // get the inside and outside volvars
278 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
279 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
280
281 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
282 const Scalar outsideWeight = 1.0 - insideWeight;
283 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
284 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
285
286 // diffusive term
287 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
288 const auto& fluxCache = elemFluxVarsCache[scvf];
289 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
290 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
291 const Scalar massOrMolarDensity = 0.5*(rhoInside + rhoOutside);
292
293 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
294 {
295 // diffusive term
296 Scalar diffDeriv = 0.0;
297 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
298 {
299 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
300 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx);
301 }
302 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
303 {
304 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)
305 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
306 }
307 else
308 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
309
310 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
311 if (!scvf.boundary())
312 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
313 }
314 }
315
316 template<class JacobianMatrix, class T = TypeTag>
317 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
318 addFluxDerivatives(JacobianMatrix& A,
319 const Problem& problem,
320 const Element& element,
321 const FVElementGeometry& fvGeometry,
322 const ElementVolumeVariables& curElemVolVars,
323 const ElementFluxVariablesCache& elemFluxVarsCache,
324 const SubControlVolumeFace& scvf) const
325 {
326
327 // advective term: we do the same for all tracer components
328 auto rho = [](const VolumeVariables& volVars)
329 { return useMoles ? volVars.molarDensity() : volVars.density(); };
330
331 // the volume flux
332 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
333
334 // the upwind weight
335 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
336
337 // get the inside and outside volvars
338 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
339 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
340
341 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
342 const auto outsideWeight = 1.0 - insideWeight;
343 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
344 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
345
346 // diffusive term
347 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
349 const auto ti = DiffusionType::calculateTransmissibilities(problem,
350 element,
351 fvGeometry,
352 curElemVolVars,
353 scvf,
354 elemFluxVarsCache[scvf],
355 phaseIdx);
356 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
357 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
358
359 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
360 {
361 for (const auto& scv : scvs(fvGeometry))
362 {
363 // diffusive term
364 auto diffDeriv = 0.0;
365 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
366 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
367 : ti[compIdx][scv.indexInElement()];
368 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
369 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
370 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
371 else
372 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
373 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
374 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
375 }
376
377 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
378 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
379 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
380 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
381 }
382 }
383
384 template<class PartialDerivativeMatrices>
385 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
386 const Problem& problem,
387 const Element& element,
388 const FVElementGeometry& fvGeometry,
389 const ElementVolumeVariables& curElemVolVars,
390 const ElementFluxVariablesCache& elemFluxVarsCache,
391 const SubControlVolumeFace& scvf) const
392 {
393 // do the same as for inner facets
394 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry,
395 curElemVolVars, elemFluxVarsCache, scvf);
396 }
397
398 template<class PartialDerivativeMatrices>
399 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
400 const Problem& problem,
401 const Element& element,
402 const FVElementGeometry& fvGeometry,
403 const ElementVolumeVariables& curElemVolVars,
404 const ElementFluxVariablesCache& elemFluxVarsCache,
405 const SubControlVolumeFace& scvf) const
406 {
407 // TODO maybe forward to the problem?
408 }
409
410private:
411 template <class T = ModelTraits>
412 [[deprecated("All compositional models must specifiy if dispersion is enabled."
413 "Please add enableCompositionalDispersion to the ModelTraits in your model header.")]]
414 void enableCompositionalDispersionMissing_() const {}
415};
416
417} // end namespace Dumux
418
419#endif
A helper to deduce a vector with the same size as numbers of equations.
Helpers for deprecation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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:55
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 Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
constexpr CCTpfa cctpfa
Definition: method.hh:137
constexpr Box box
Definition: method.hh:139
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
std::string porosity() noexcept
I/O name of porosity.
Definition: name.hh:139
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Element-wise calculation of the local residual for problems using fully implicit tracer model.
Definition: porousmediumflow/tracer/localresidual.hh:48
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:85
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:399
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:256
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:244
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:212
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:130
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:318
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:385
Declares all properties used in Dumux.