26#ifndef DUMUX_TRACER_LOCAL_RESIDUAL_HH
27#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
29#include <dune/common/exceptions.hh>
45template<
class TypeTag>
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolume =
typename GridGeometry::SubControlVolume;
54 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
60 using Element =
typename GridView::template Codim<0>::Entity;
66 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
67 static constexpr int phaseIdx = 0;
70 using ParentType::ParentType;
84 const SubControlVolume& scv,
85 const VolumeVariables& volVars)
const
87 NumEqVector storage(0.0);
93 const Scalar
saturation = max(1e-8, volVars.saturation(phaseIdx));
98 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
99 storage[compIdx] += volVars.porosity()
100 * volVars.molarDensity(phaseIdx)
101 * volVars.moleFraction(phaseIdx, compIdx)
107 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
108 storage[compIdx] += volVars.porosity()
109 * volVars.density(phaseIdx)
110 * volVars.massFraction(phaseIdx, compIdx)
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& elemVolVars,
132 const SubControlVolumeFace& scvf,
133 const ElementFluxVariablesCache& elemFluxVarsCache)
const
135 FluxVariables fluxVars;
136 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
138 NumEqVector flux(0.0);
139 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
141 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
145 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
148 auto upwindTerm = [compIdx](
const auto& volVars)
149 {
return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
152 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
155 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
157 flux[compIdx] += diffusiveFluxes[compIdx];
159 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
165 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
168 auto upwindTerm = [compIdx](
const auto& volVars)
169 {
return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
172 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
175 flux[compIdx] += diffusiveFluxes[compIdx];
177 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
179 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
196 template<
class PartialDerivativeMatrix>
198 const Problem& problem,
199 const Element& element,
200 const FVElementGeometry& fvGeometry,
201 const VolumeVariables& curVolVars,
202 const SubControlVolume& scv)
const
208 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
210 const auto porosity = curVolVars.porosity();
211 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
214 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
215 partialDerivatives[compIdx][compIdx] += d_storage;
228 template<
class PartialDerivativeMatrix>
230 const Problem& problem,
231 const Element& element,
232 const FVElementGeometry& fvGeometry,
233 const VolumeVariables& curVolVars,
234 const SubControlVolume& scv)
const
239 template<
class PartialDerivativeMatrices,
class T = TypeTag>
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
250 DUNE_THROW(Dune::NotImplemented,
"Analytic flux differentiation only implemented for tpfa");
253 auto rho = [](
const VolumeVariables& volVars)
254 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
257 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
260 static const Scalar upwindWeight = getParam<Scalar>(
"Flux.UpwindWeight");
263 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
264 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
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);
278 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
281 Scalar diffDeriv = 0.0;
284 diffDeriv = useMoles ?
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
290 :
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
293 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
295 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
296 if (!scvf.boundary())
297 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
301 template<
class JacobianMatrix,
class T = TypeTag>
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
313 auto rho = [](
const VolumeVariables& volVars)
314 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
317 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
320 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
323 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
324 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
332 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
334 const auto ti = DiffusionType::calculateTransmissibilities(problem,
339 elemFluxVarsCache[scvf],
341 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
342 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
344 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
346 for (
const auto& scv : scvs(fvGeometry))
349 auto diffDeriv = 0.0;
351 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
352 : ti[compIdx][scv.indexInElement()];
354 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
355 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
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;
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;
369 template<
class PartialDerivativeMatrices>
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
380 curElemVolVars, elemFluxVarsCache, scvf);
383 template<
class PartialDerivativeMatrices>
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
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
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.