3.1-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
33
34namespace Dumux {
35
41template<class TypeTag>
42class TracerLocalResidual: public GetPropType<TypeTag, Properties::BaseLocalResidual>
43{
47 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
48 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
49 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
52 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
54 using Element = typename GridView::template Codim<0>::Entity;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
58
59 static constexpr int numComponents = GetPropType<TypeTag, Properties::ModelTraits>::numFluidComponents();
60 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
61 static constexpr int phaseIdx = 0;
62
63public:
64 using ParentType::ParentType;
65
77 NumEqVector computeStorage(const Problem& problem,
78 const SubControlVolume& scv,
79 const VolumeVariables& volVars) const
80 {
81 NumEqVector storage(0.0);
82
83 // formulation with mole balances
84 if (useMoles)
85 {
86 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
87 storage[compIdx] += volVars.porosity()
88 * volVars.molarDensity(phaseIdx)
89 * volVars.moleFraction(phaseIdx, compIdx)
90 * volVars.saturation(phaseIdx);
91 }
92 // formulation with mass balances
93 else
94 {
95 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
96 storage[compIdx] += volVars.porosity()
97 * volVars.density(phaseIdx)
98 * volVars.massFraction(phaseIdx, compIdx)
99 * volVars.saturation(phaseIdx);
100
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 // formulation with mole balances
132 if (useMoles)
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 // formulation with mass balances
152 else
153 {
154 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
155 {
156 // the physical quantities for which we perform upwinding
157 auto upwindTerm = [compIdx](const auto& volVars)
158 { return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
159
160 // advective fluxes
161 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
162 // diffusive fluxes
163 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
164 flux[compIdx] += diffusiveFluxes[compIdx];
165 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
166 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
167 else
168 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
169 }
170 }
171
172 return flux;
173 }
174
185 template<class PartialDerivativeMatrix>
186 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
187 const Problem& problem,
188 const Element& element,
189 const FVElementGeometry& fvGeometry,
190 const VolumeVariables& curVolVars,
191 const SubControlVolume& scv) const
192 {
193 const auto porosity = curVolVars.porosity();
194 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
195 const auto d_storage = scv.volume()*porosity*rho/this->timeLoop().timeStepSize();
196
197 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
198 partialDerivatives[compIdx][compIdx] += d_storage;
199 }
200
211 template<class PartialDerivativeMatrix>
212 void addSourceDerivatives(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 // TODO maybe forward to the problem? -> necessary for reaction terms
220 }
221
222 template<class PartialDerivativeMatrices, class T = TypeTag>
223 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod != DiscretizationMethod::box, void>
224 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
225 const Problem& problem,
226 const Element& element,
227 const FVElementGeometry& fvGeometry,
228 const ElementVolumeVariables& curElemVolVars,
229 const ElementFluxVariablesCache& elemFluxVarsCache,
230 const SubControlVolumeFace& scvf) const
231 {
232 // advective term: we do the same for all tracer components
233 auto rho = [](const VolumeVariables& volVars)
234 { return useMoles ? volVars.molarDensity() : volVars.density(); };
235
236 // the volume flux
237 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
238
239 // the upwind weight
240 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
241
242 // get the inside and outside volvars
243 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
244 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
245
246 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
247 const Scalar outsideWeight = 1.0 - insideWeight;
248 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
249 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
250
251 // diffusive term
252 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
253 const auto& fluxCache = elemFluxVarsCache[scvf];
254 const Scalar rhoInside = massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
255 const Scalar rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
256 const Scalar massOrMolarDensity = 0.5*(rhoInside + rhoOutside);
257
258 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
259 {
260 // diffusive term
261 Scalar diffDeriv = 0.0;
262 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
263 {
264 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
265 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx);
266 }
267 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
268 {
269 diffDeriv = useMoles ? massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)
270 : massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
271 }
272 else
273 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
274
275 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
276 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
277 }
278 }
279
280 template<class JacobianMatrix, class T = TypeTag>
281 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::box, void>
282 addFluxDerivatives(JacobianMatrix& A,
283 const Problem& problem,
284 const Element& element,
285 const FVElementGeometry& fvGeometry,
286 const ElementVolumeVariables& curElemVolVars,
287 const ElementFluxVariablesCache& elemFluxVarsCache,
288 const SubControlVolumeFace& scvf) const
289 {
290
291 // advective term: we do the same for all tracer components
292 auto rho = [](const VolumeVariables& volVars)
293 { return useMoles ? volVars.molarDensity() : volVars.density(); };
294
295 // the volume flux
296 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
297
298 // the upwind weight
299 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
300
301 // get the inside and outside volvars
302 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
303 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
304
305 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
306 const auto outsideWeight = 1.0 - insideWeight;
307 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
308 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
309
310 // diffusive term
311 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
313 const auto ti = DiffusionType::calculateTransmissibilities(problem,
314 element,
315 fvGeometry,
316 curElemVolVars,
317 scvf,
318 elemFluxVarsCache[scvf],
319 phaseIdx);
320 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
321 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
322
323 for (int compIdx = 0; compIdx < numComponents; ++compIdx)
324 {
325 for (const auto& scv : scvs(fvGeometry))
326 {
327 // diffusive term
328 auto diffDeriv = 0.0;
329 if (referenceSystemFormulation == ReferenceSystemFormulation::massAveraged)
330 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
331 : ti[compIdx][scv.indexInElement()];
332 else if (referenceSystemFormulation == ReferenceSystemFormulation::molarAveraged)
333 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
334 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
335 else
336 DUNE_THROW(Dune::NotImplemented, "other reference systems than mass and molar averaged are not implemented");
337 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
338 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
339 }
340
341 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
342 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
343 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
344 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
345 }
346 }
347
348 template<class PartialDerivativeMatrices>
349 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
350 const Problem& problem,
351 const Element& element,
352 const FVElementGeometry& fvGeometry,
353 const ElementVolumeVariables& curElemVolVars,
354 const ElementFluxVariablesCache& elemFluxVarsCache,
355 const SubControlVolumeFace& scvf) const
356 {
357 // do the same as for inner facets
358 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry,
359 curElemVolVars, elemFluxVarsCache, scvf);
360 }
361
362 template<class PartialDerivativeMatrices>
363 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
364 const Problem& problem,
365 const Element& element,
366 const FVElementGeometry& fvGeometry,
367 const ElementVolumeVariables& curElemVolVars,
368 const ElementFluxVariablesCache& elemFluxVarsCache,
369 const SubControlVolumeFace& scvf) const
370 {
371 // TODO maybe forward to the problem?
372 }
373};
374
375} // end namespace Dumux
376
377#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
make the local view function available whenever we use the grid geometry
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 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:43
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:77
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:282
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:363
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:212
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:186
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 !=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:224
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:349
Declares all properties used in Dumux.