3.4
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
37
38namespace Dumux {
39
45template<class TypeTag>
46class TracerLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
47{
52 using FVElementGeometry = typename GridGeometry::LocalView;
53 using SubControlVolume = typename GridGeometry::SubControlVolume;
54 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
55 using Extrusion = Extrusion_t<GridGeometry>;
58 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
60 using Element = typename GridView::template Codim<0>::Entity;
61 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
64
65 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
66 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
67 static constexpr int phaseIdx = 0;
68
69public:
70 using ParentType::ParentType;
71
83 NumEqVector computeStorage(const Problem& problem,
84 const SubControlVolume& scv,
85 const VolumeVariables& volVars) const
86 {
87 NumEqVector storage(0.0);
88
89 // regularize saturation so we don't get singular matrices when the saturation is zero
90 // note that the fluxes will still be zero (zero effective diffusion coefficient),
91 // and we still solve the equation storage = 0 yielding the correct result
92 using std::max;
93 const Scalar saturation = max(1e-8, volVars.saturation(phaseIdx));
94
95 // formulation with mole balances
96 if (useMoles)
97 {
98 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
99 storage[compIdx] += volVars.porosity()
100 * volVars.molarDensity(phaseIdx)
101 * volVars.moleFraction(phaseIdx, compIdx)
102 * saturation;
103 }
104 // formulation with mass balances
105 else
106 {
107 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
108 storage[compIdx] += volVars.porosity()
109 * volVars.density(phaseIdx)
110 * volVars.massFraction(phaseIdx, compIdx)
111 * saturation;
112 }
113
114 return storage;
115 }
116
128 NumEqVector computeFlux(const Problem& problem,
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const SubControlVolumeFace& scvf,
133 const ElementFluxVariablesCache& elemFluxVarsCache) const
134 {
135 FluxVariables fluxVars;
136 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
137
138 NumEqVector flux(0.0);
139 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
140
141 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
142 // formulation with mole balances
143 if (useMoles)
144 {
145 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
146 {
147 // the physical quantities for which we perform upwinding
148 auto upwindTerm = [compIdx](const auto& volVars)
149 { return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
150
151 // advective fluxes
152 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
153 // diffusive fluxes
154 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
155 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
156 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
157 flux[compIdx] += diffusiveFluxes[compIdx];
158 else
159 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
160 }
161 }
162 // formulation with mass balances
163 else
164 {
165 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
166 {
167 // the physical quantities for which we perform upwinding
168 auto upwindTerm = [compIdx](const auto& volVars)
169 { return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
170
171 // advective fluxes
172 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
173 // diffusive fluxes
174 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
175 flux[compIdx] += diffusiveFluxes[compIdx];
176 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
177 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
178 else
179 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
180 }
181 }
182
183 return flux;
184 }
185
196 template<class PartialDerivativeMatrix>
197 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
198 const Problem& problem,
199 const Element& element,
200 const FVElementGeometry& fvGeometry,
201 const VolumeVariables& curVolVars,
202 const SubControlVolume& scv) const
203 {
204 // regularize saturation so we don't get singular matrices when the saturation is zero
205 // note that the fluxes will still be zero (zero effective diffusion coefficient),
206 // and we still solve the equation storage = 0 yielding the correct result
207 using std::max;
208 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
209
210 const auto porosity = curVolVars.porosity();
211 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
212 const auto d_storage = Extrusion::volume(scv)*porosity*rho*saturation/this->timeLoop().timeStepSize();
213
214 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
215 partialDerivatives[compIdx][compIdx] += d_storage;
216 }
217
228 template<class PartialDerivativeMatrix>
229 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
230 const Problem& problem,
231 const Element& element,
232 const FVElementGeometry& fvGeometry,
233 const VolumeVariables& curVolVars,
234 const SubControlVolume& scv) const
235 {
236 // TODO maybe forward to the problem? -> necessary for reaction terms
237 }
238
239 template<class PartialDerivativeMatrices, class T = TypeTag>
240 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethod::box, void>
241 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
242 const Problem& problem,
243 const Element& element,
244 const FVElementGeometry& fvGeometry,
245 const ElementVolumeVariables& curElemVolVars,
246 const ElementFluxVariablesCache& elemFluxVarsCache,
247 const SubControlVolumeFace& scvf) const
248 {
249 if constexpr (FVElementGeometry::GridGeometry::discMethod != DiscretizationMethod::cctpfa)
250 DUNE_THROW(Dune::NotImplemented, "Analytic flux differentiation only implemented for tpfa");
251
252 // advective term: we do the same for all tracer components
253 auto rho = [](const VolumeVariables& volVars)
254 { return useMoles ? volVars.molarDensity() : volVars.density(); };
255
256 // the volume flux
257 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
258
259 // the upwind weight
260 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
261
262 // get the inside and outside volvars
263 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
264 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
265
266 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
267 const Scalar outsideWeight = 1.0 - insideWeight;
268 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
269 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
270
271 // diffusive term
272 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
273 const auto& fluxCache = elemFluxVarsCache[scvf];
274 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
275 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
276 const Scalar massOrMolarDensity = 0.5*(rhoInside + rhoOutside);
277
278 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
279 {
280 // diffusive term
281 Scalar diffDeriv = 0.0;
282 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
283 {
284 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
285 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx);
286 }
287 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
288 {
289 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)
290 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
291 }
292 else
293 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
294
295 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
296 if (!scvf.boundary())
297 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
298 }
299 }
300
301 template<class JacobianMatrix, class T = TypeTag>
302 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::box, void>
303 addFluxDerivatives(JacobianMatrix& A,
304 const Problem& problem,
305 const Element& element,
306 const FVElementGeometry& fvGeometry,
307 const ElementVolumeVariables& curElemVolVars,
308 const ElementFluxVariablesCache& elemFluxVarsCache,
309 const SubControlVolumeFace& scvf) const
310 {
311
312 // advective term: we do the same for all tracer components
313 auto rho = [](const VolumeVariables& volVars)
314 { return useMoles ? volVars.molarDensity() : volVars.density(); };
315
316 // the volume flux
317 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
318
319 // the upwind weight
320 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
321
322 // get the inside and outside volvars
323 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
324 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
325
326 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
327 const auto outsideWeight = 1.0 - insideWeight;
328 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
329 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
330
331 // diffusive term
332 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
334 const auto ti = DiffusionType::calculateTransmissibilities(problem,
335 element,
336 fvGeometry,
337 curElemVolVars,
338 scvf,
339 elemFluxVarsCache[scvf],
340 phaseIdx);
341 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
342 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
343
344 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
345 {
346 for (const auto& scv : scvs(fvGeometry))
347 {
348 // diffusive term
349 auto diffDeriv = 0.0;
350 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
351 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
352 : ti[compIdx][scv.indexInElement()];
353 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
354 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
355 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
356 else
357 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
358 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
359 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
360 }
361
362 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
363 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
364 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
365 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
366 }
367 }
368
369 template<class PartialDerivativeMatrices>
370 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
371 const Problem& problem,
372 const Element& element,
373 const FVElementGeometry& fvGeometry,
374 const ElementVolumeVariables& curElemVolVars,
375 const ElementFluxVariablesCache& elemFluxVarsCache,
376 const SubControlVolumeFace& scvf) const
377 {
378 // do the same as for inner facets
379 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry,
380 curElemVolVars, elemFluxVarsCache, scvf);
381 }
382
383 template<class PartialDerivativeMatrices>
384 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
385 const Problem& problem,
386 const Element& element,
387 const FVElementGeometry& fvGeometry,
388 const ElementVolumeVariables& curElemVolVars,
389 const ElementFluxVariablesCache& elemFluxVarsCache,
390 const SubControlVolumeFace& scvf) const
391 {
392 // TODO maybe forward to the problem?
393 }
394};
395
396} // end namespace Dumux
397
398#endif
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 available discretization methods in Dumux.
Helper classes to compute the integration elements.
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
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:47
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:83
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::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:303
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:384
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:229
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:197
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:128
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod !=DiscretizationMethod::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:241
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:370
Declares all properties used in Dumux.