3.6-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
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;
65
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
68 static constexpr int phaseIdx = 0;
69
70public:
71 using ParentType::ParentType;
72
84 NumEqVector computeStorage(const Problem& problem,
85 const SubControlVolume& scv,
86 const VolumeVariables& volVars) const
87 {
88 NumEqVector storage(0.0);
89
90 // regularize saturation so we don't get singular matrices when the saturation is zero
91 // note that the fluxes will still be zero (zero effective diffusion coefficient),
92 // and we still solve the equation storage = 0 yielding the correct result
93 using std::max;
94 const Scalar saturation = max(1e-8, volVars.saturation(phaseIdx));
95
96 // formulation with mole balances
97 if (useMoles)
98 {
99 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
100 storage[compIdx] += volVars.porosity()
101 * volVars.molarDensity(phaseIdx)
102 * volVars.moleFraction(phaseIdx, compIdx)
103 * saturation;
104 }
105 // formulation with mass balances
106 else
107 {
108 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
109 storage[compIdx] += volVars.porosity()
110 * volVars.density(phaseIdx)
111 * volVars.massFraction(phaseIdx, compIdx)
112 * saturation;
113 }
114
115 return storage;
116 }
117
129 NumEqVector computeFlux(const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache) const
135 {
136 FluxVariables fluxVars;
137 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
138
139 NumEqVector flux(0.0);
140 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
141
142 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
143
144 if (useMoles) // formulation with mole balances
145 {
146 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
147 {
148 // the physical quantities for which we perform upwinding
149 auto upwindTerm = [compIdx](const auto& volVars)
150 { return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
151
152 // advective fluxes
153 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
154 // diffusive fluxes
155 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
156 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
157 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
158 flux[compIdx] += diffusiveFluxes[compIdx];
159 else
160 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
161 }
162 }
163 else // formulation with mass balances
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 if constexpr (ModelTraits::enableCompositionalDispersion())
184 {
185 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
186 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
187 {
188 flux[compIdx] += dispersionFluxes[compIdx];
189 }
190 }
191
192 return flux;
193 }
194
205 template<class PartialDerivativeMatrix>
206 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
207 const Problem& problem,
208 const Element& element,
209 const FVElementGeometry& fvGeometry,
210 const VolumeVariables& curVolVars,
211 const SubControlVolume& scv) const
212 {
213 // regularize saturation so we don't get singular matrices when the saturation is zero
214 // note that the fluxes will still be zero (zero effective diffusion coefficient),
215 // and we still solve the equation storage = 0 yielding the correct result
216 using std::max;
217 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
218
219 const auto porosity = curVolVars.porosity();
220 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
221 const auto d_storage = Extrusion::volume(fvGeometry, scv)*porosity*rho*saturation/this->timeLoop().timeStepSize();
222
223 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
224 partialDerivatives[compIdx][compIdx] += d_storage;
225 }
226
237 template<class PartialDerivativeMatrix>
238 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
239 const Problem& problem,
240 const Element& element,
241 const FVElementGeometry& fvGeometry,
242 const VolumeVariables& curVolVars,
243 const SubControlVolume& scv) const
244 {
245 // TODO maybe forward to the problem? -> necessary for reaction terms
246 }
247
248 template<class PartialDerivativeMatrices, class T = TypeTag>
249 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethods::box, void>
250 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
251 const Problem& problem,
252 const Element& element,
253 const FVElementGeometry& fvGeometry,
254 const ElementVolumeVariables& curElemVolVars,
255 const ElementFluxVariablesCache& elemFluxVarsCache,
256 const SubControlVolumeFace& scvf) const
257 {
258 if constexpr (FVElementGeometry::GridGeometry::discMethod != DiscretizationMethods::cctpfa)
259 DUNE_THROW(Dune::NotImplemented, "Analytic flux differentiation only implemented for tpfa");
260
261 // advective term: we do the same for all tracer components
262 auto rho = [](const VolumeVariables& volVars)
263 { return useMoles ? volVars.molarDensity() : volVars.density(); };
264
265 // the volume flux
266 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
267
268 // the upwind weight
269 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
270
271 // get the inside and outside volvars
272 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
273 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
274
275 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
276 const Scalar outsideWeight = 1.0 - insideWeight;
277 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
278 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
279
280 // diffusive term
281 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
282 const auto& fluxCache = elemFluxVarsCache[scvf];
283 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
284 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
285 const Scalar massOrMolarDensity = 0.5*(rhoInside + rhoOutside);
286
287 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
288 {
289 // diffusive term
290 Scalar diffDeriv = 0.0;
291 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
292 {
293 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
294 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx);
295 }
296 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
297 {
298 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)
299 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
300 }
301 else
302 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
303
304 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
305 if (!scvf.boundary())
306 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
307 }
308 }
309
310 template<class JacobianMatrix, class T = TypeTag>
311 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
312 addFluxDerivatives(JacobianMatrix& A,
313 const Problem& problem,
314 const Element& element,
315 const FVElementGeometry& fvGeometry,
316 const ElementVolumeVariables& curElemVolVars,
317 const ElementFluxVariablesCache& elemFluxVarsCache,
318 const SubControlVolumeFace& scvf) const
319 {
320
321 // advective term: we do the same for all tracer components
322 auto rho = [](const VolumeVariables& volVars)
323 { return useMoles ? volVars.molarDensity() : volVars.density(); };
324
325 // the volume flux
326 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
327
328 // the upwind weight
329 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
330
331 // get the inside and outside volvars
332 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
333 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
334
335 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
336 const auto outsideWeight = 1.0 - insideWeight;
337 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
338 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
339
340 // diffusive term
341 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
343 const auto ti = DiffusionType::calculateTransmissibilities(problem,
344 element,
345 fvGeometry,
346 curElemVolVars,
347 scvf,
348 elemFluxVarsCache[scvf],
349 phaseIdx);
350 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
351 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
352
353 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
354 {
355 for (const auto& scv : scvs(fvGeometry))
356 {
357 // diffusive term
358 auto diffDeriv = 0.0;
359 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
360 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
361 : ti[compIdx][scv.indexInElement()];
362 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
363 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
364 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
365 else
366 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
367 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
368 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
369 }
370
371 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
372 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
373 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
374 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
375 }
376 }
377
378 template<class PartialDerivativeMatrices>
379 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
380 const Problem& problem,
381 const Element& element,
382 const FVElementGeometry& fvGeometry,
383 const ElementVolumeVariables& curElemVolVars,
384 const ElementFluxVariablesCache& elemFluxVarsCache,
385 const SubControlVolumeFace& scvf) const
386 {
387 // do the same as for inner facets
388 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry,
389 curElemVolVars, elemFluxVarsCache, scvf);
390 }
391
392 template<class PartialDerivativeMatrices>
393 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
394 const Problem& problem,
395 const Element& element,
396 const FVElementGeometry& fvGeometry,
397 const ElementVolumeVariables& curElemVolVars,
398 const ElementFluxVariablesCache& elemFluxVarsCache,
399 const SubControlVolumeFace& scvf) const
400 {
401 // TODO maybe forward to the problem?
402 }
403};
404
405} // end namespace Dumux
406
407#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.
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...
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr CCTpfa cctpfa
Definition: method.hh:134
constexpr Box box
Definition: method.hh:136
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
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:84
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:393
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:250
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:238
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:206
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:129
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:312
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:379
Declares all properties used in Dumux.