26#ifndef DUMUX_TRACER_LOCAL_RESIDUAL_HH
27#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
29#include <dune/common/exceptions.hh>
44template<
class TypeTag>
51 using FVElementGeometry =
typename GridGeometry::LocalView;
52 using SubControlVolume =
typename GridGeometry::SubControlVolume;
53 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
59 using Element =
typename GridView::template Codim<0>::Entity;
65 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
66 static constexpr int phaseIdx = 0;
69 using ParentType::ParentType;
83 const SubControlVolume& scv,
84 const VolumeVariables& volVars)
const
86 NumEqVector storage(0.0);
92 const Scalar
saturation = max(1e-8, volVars.saturation(phaseIdx));
97 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
98 storage[compIdx] += volVars.porosity()
99 * volVars.molarDensity(phaseIdx)
100 * volVars.moleFraction(phaseIdx, compIdx)
106 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
107 storage[compIdx] += volVars.porosity()
108 * volVars.density(phaseIdx)
109 * volVars.massFraction(phaseIdx, compIdx)
128 const Element& element,
129 const FVElementGeometry& fvGeometry,
130 const ElementVolumeVariables& elemVolVars,
131 const SubControlVolumeFace& scvf,
132 const ElementFluxVariablesCache& elemFluxVarsCache)
const
134 FluxVariables fluxVars;
135 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
137 NumEqVector flux(0.0);
138 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
140 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
144 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
147 auto upwindTerm = [compIdx](
const auto& volVars)
148 {
return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
151 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
154 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
156 flux[compIdx] += diffusiveFluxes[compIdx];
158 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
164 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
167 auto upwindTerm = [compIdx](
const auto& volVars)
168 {
return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
171 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
174 flux[compIdx] += diffusiveFluxes[compIdx];
176 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
178 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
195 template<
class PartialDerivativeMatrix>
197 const Problem& problem,
198 const Element& element,
199 const FVElementGeometry& fvGeometry,
200 const VolumeVariables& curVolVars,
201 const SubControlVolume& scv)
const
207 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
209 const auto porosity = curVolVars.porosity();
210 const auto rho = useMoles ? curVolVars.molarDensity() : curVolVars.density();
211 const auto d_storage = Extrusion::volume(scv)*
porosity*rho*
saturation/this->timeLoop().timeStepSize();
213 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
214 partialDerivatives[compIdx][compIdx] += d_storage;
227 template<
class PartialDerivativeMatrix>
229 const Problem& problem,
230 const Element& element,
231 const FVElementGeometry& fvGeometry,
232 const VolumeVariables& curVolVars,
233 const SubControlVolume& scv)
const
238 template<
class PartialDerivativeMatrices,
class T = TypeTag>
241 const Problem& problem,
242 const Element& element,
243 const FVElementGeometry& fvGeometry,
244 const ElementVolumeVariables& curElemVolVars,
245 const ElementFluxVariablesCache& elemFluxVarsCache,
246 const SubControlVolumeFace& scvf)
const
249 DUNE_THROW(Dune::NotImplemented,
"Analytic flux differentiation only implemented for tpfa");
252 auto rho = [](
const VolumeVariables& volVars)
253 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
256 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
259 static const Scalar upwindWeight = getParam<Scalar>(
"Flux.UpwindWeight");
262 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
263 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
265 const Scalar insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
266 const Scalar outsideWeight = 1.0 - insideWeight;
267 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
268 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
271 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
272 const auto& fluxCache = elemFluxVarsCache[scvf];
273 const Scalar rhoInside =
massOrMolarDensity(insideVolVars, referenceSystemFormulation, phaseIdx);
274 const Scalar rhoOutside =
massOrMolarDensity(outsideVolVars, referenceSystemFormulation, phaseIdx);
277 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
280 Scalar diffDeriv = 0.0;
283 diffDeriv = useMoles ?
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
289 :
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
292 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
294 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
295 if (!scvf.boundary())
296 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
300 template<
class JacobianMatrix,
class T = TypeTag>
303 const Problem& problem,
304 const Element& element,
305 const FVElementGeometry& fvGeometry,
306 const ElementVolumeVariables& curElemVolVars,
307 const ElementFluxVariablesCache& elemFluxVarsCache,
308 const SubControlVolumeFace& scvf)
const
312 auto rho = [](
const VolumeVariables& volVars)
313 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
316 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
319 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
322 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
323 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
325 const auto insideWeight = std::signbit(volFlux) ? (1.0 - upwindWeight) : upwindWeight;
326 const auto outsideWeight = 1.0 - insideWeight;
327 const auto advDerivII = volFlux*rho(insideVolVars)*insideWeight;
328 const auto advDerivIJ = volFlux*rho(outsideVolVars)*outsideWeight;
331 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
333 const auto ti = DiffusionType::calculateTransmissibilities(problem,
338 elemFluxVarsCache[scvf],
340 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
341 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
343 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
345 for (
const auto& scv : scvs(fvGeometry))
348 auto diffDeriv = 0.0;
350 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
351 : ti[compIdx][scv.indexInElement()];
353 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
354 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
356 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
357 A[insideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] += diffDeriv;
358 A[outsideScv.dofIndex()][scv.dofIndex()][compIdx][compIdx] -= diffDeriv;
361 A[insideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] += advDerivII;
362 A[insideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] += advDerivIJ;
363 A[outsideScv.dofIndex()][outsideScv.dofIndex()][compIdx][compIdx] -= advDerivII;
364 A[outsideScv.dofIndex()][insideScv.dofIndex()][compIdx][compIdx] -= advDerivIJ;
368 template<
class PartialDerivativeMatrices>
370 const Problem& problem,
371 const Element& element,
372 const FVElementGeometry& fvGeometry,
373 const ElementVolumeVariables& curElemVolVars,
374 const ElementFluxVariablesCache& elemFluxVarsCache,
375 const SubControlVolumeFace& scvf)
const
379 curElemVolVars, elemFluxVarsCache, scvf);
382 template<
class PartialDerivativeMatrices>
384 const Problem& problem,
385 const Element& element,
386 const FVElementGeometry& fvGeometry,
387 const ElementVolumeVariables& curElemVolVars,
388 const ElementFluxVariablesCache& elemFluxVarsCache,
389 const SubControlVolumeFace& scvf)
const
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 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 (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:46
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:82
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:302
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:383
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:228
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:196
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:127
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:240
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:369
Declares all properties used in Dumux.