3.2-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
35
36namespace Dumux {
37
43template<class TypeTag>
44class TracerLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
45{
49 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
50 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
54 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
56 using Element = typename GridView::template Codim<0>::Entity;
57 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
60
61 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
62 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
63 static constexpr int phaseIdx = 0;
64
65public:
66 using ParentType::ParentType;
67
79 NumEqVector computeStorage(const Problem& problem,
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars) const
82 {
83 NumEqVector storage(0.0);
84
85 // regularize saturation so we don't get singular matrices when the saturation is zero
86 // note that the fluxes will still be zero (zero effective diffusion coefficient),
87 // and we still solve the equation storage = 0 yielding the correct result
88 using std::max;
89 const Scalar saturation = max(1e-8, volVars.saturation(phaseIdx));
90
91 // formulation with mole balances
92 if (useMoles)
93 {
94 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
95 storage[compIdx] += volVars.porosity()
96 * volVars.molarDensity(phaseIdx)
97 * volVars.moleFraction(phaseIdx, compIdx)
98 * saturation;
99 }
100 // formulation with mass balances
101 else
102 {
103 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
104 storage[compIdx] += volVars.porosity()
105 * volVars.density(phaseIdx)
106 * volVars.massFraction(phaseIdx, compIdx)
107 * saturation;
108 }
109
110 return storage;
111 }
112
124 NumEqVector computeFlux(const Problem& problem,
125 const Element& element,
126 const FVElementGeometry& fvGeometry,
127 const ElementVolumeVariables& elemVolVars,
128 const SubControlVolumeFace& scvf,
129 const ElementFluxVariablesCache& elemFluxVarsCache) const
130 {
131 FluxVariables fluxVars;
132 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
133
134 NumEqVector flux(0.0);
135 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
136
137 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
138 // formulation with mole balances
139 if (useMoles)
140 {
141 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
142 {
143 // the physical quantities for which we perform upwinding
144 auto upwindTerm = [compIdx](const auto& volVars)
145 { return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
146
147 // advective fluxes
148 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
149 // diffusive fluxes
150 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
151 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
152 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
153 flux[compIdx] += diffusiveFluxes[compIdx];
154 else
155 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
156 }
157 }
158 // formulation with mass balances
159 else
160 {
161 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
162 {
163 // the physical quantities for which we perform upwinding
164 auto upwindTerm = [compIdx](const auto& volVars)
165 { return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
166
167 // advective fluxes
168 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
169 // diffusive fluxes
170 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
171 flux[compIdx] += diffusiveFluxes[compIdx];
172 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
173 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
174 else
175 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
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 = scv.volume()*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 != DiscretizationMethod::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 != DiscretizationMethod::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 == DiscretizationMethod::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
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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:45
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:79
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:299
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
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:124
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:237
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
Declares all properties used in Dumux.