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 int numComponents = ModelTraits::numFluidComponents();
67 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
68 static constexpr int phaseIdx = 0;
71 using ParentType::ParentType;
85 const SubControlVolume& scv,
86 const VolumeVariables& volVars)
const
88 NumEqVector storage(0.0);
94 const Scalar
saturation = max(1e-8, volVars.saturation(phaseIdx));
99 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
100 storage[compIdx] += volVars.porosity()
101 * volVars.molarDensity(phaseIdx)
102 * volVars.moleFraction(phaseIdx, compIdx)
108 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
109 storage[compIdx] += volVars.porosity()
110 * volVars.density(phaseIdx)
111 * volVars.massFraction(phaseIdx, compIdx)
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache)
const
136 FluxVariables fluxVars;
137 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
139 NumEqVector flux(0.0);
140 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
142 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
146 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
149 auto upwindTerm = [compIdx](
const auto& volVars)
150 {
return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
153 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
156 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
158 flux[compIdx] += diffusiveFluxes[compIdx];
160 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");
183 if constexpr (ModelTraits::enableCompositionalDispersion())
185 const auto dispersionFluxes = fluxVars.compositionalDispersionFlux(phaseIdx);
186 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
188 flux[compIdx] += dispersionFluxes[compIdx];
205 template<
class PartialDerivativeMatrix>
207 const Problem& problem,
208 const Element& element,
209 const FVElementGeometry& fvGeometry,
210 const VolumeVariables& curVolVars,
211 const SubControlVolume& scv)
const
217 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
219 const auto porosity = curVolVars.porosity();
220 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
223 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
224 partialDerivatives[compIdx][compIdx] += d_storage;
237 template<
class PartialDerivativeMatrix>
239 const Problem& problem,
240 const Element& element,
241 const FVElementGeometry& fvGeometry,
242 const VolumeVariables& curVolVars,
243 const SubControlVolume& scv)
const
248 template<
class PartialDerivativeMatrices,
class T = TypeTag>
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
259 DUNE_THROW(Dune::NotImplemented,
"Analytic flux differentiation only implemented for tpfa");
262 auto rho = [](
const VolumeVariables& volVars)
263 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
266 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
269 static const Scalar upwindWeight = getParam<Scalar>(
"Flux.UpwindWeight");
272 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
273 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
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);
287 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
290 Scalar diffDeriv = 0.0;
293 diffDeriv = useMoles ?
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
299 :
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
302 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
304 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
305 if (!scvf.boundary())
306 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
310 template<
class JacobianMatrix,
class T = TypeTag>
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
322 auto rho = [](
const VolumeVariables& volVars)
323 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
326 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
329 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
332 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
333 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
341 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
343 const auto ti = DiffusionType::calculateTransmissibilities(problem,
348 elemFluxVarsCache[scvf],
350 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
351 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
353 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
355 for (
const auto& scv : scvs(fvGeometry))
358 auto diffDeriv = 0.0;
360 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
361 : ti[compIdx][scv.indexInElement()];
363 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
364 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
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;
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;
378 template<
class PartialDerivativeMatrices>
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
389 curElemVolVars, elemFluxVarsCache, scvf);
392 template<
class PartialDerivativeMatrices>
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
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.