26#ifndef DUMUX_TRACER_LOCAL_RESIDUAL_HH
27#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
29#include <dune/common/exceptions.hh>
43template<
class TypeTag>
50 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
51 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 using Element =
typename GridView::template Codim<0>::Entity;
62 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
63 static constexpr int phaseIdx = 0;
66 using ParentType::ParentType;
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars)
const
83 NumEqVector storage(0.0);
89 const Scalar
saturation = max(1e-8, volVars.saturation(phaseIdx));
94 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
95 storage[compIdx] += volVars.porosity()
96 * volVars.molarDensity(phaseIdx)
97 * volVars.moleFraction(phaseIdx, compIdx)
103 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
104 storage[compIdx] += volVars.porosity()
105 * volVars.density(phaseIdx)
106 * volVars.massFraction(phaseIdx, compIdx)
125 const Element& element,
126 const FVElementGeometry& fvGeometry,
127 const ElementVolumeVariables& elemVolVars,
128 const SubControlVolumeFace& scvf,
129 const ElementFluxVariablesCache& elemFluxVarsCache)
const
131 FluxVariables fluxVars;
132 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
134 NumEqVector flux(0.0);
135 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
137 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
141 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
144 auto upwindTerm = [compIdx](
const auto& volVars)
145 {
return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
148 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
151 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
153 flux[compIdx] += diffusiveFluxes[compIdx];
155 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
161 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
164 auto upwindTerm = [compIdx](
const auto& volVars)
165 {
return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
168 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
171 flux[compIdx] += diffusiveFluxes[compIdx];
173 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
175 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
192 template<
class PartialDerivativeMatrix>
194 const Problem& problem,
195 const Element& element,
196 const FVElementGeometry& fvGeometry,
197 const VolumeVariables& curVolVars,
198 const SubControlVolume& scv)
const
204 const auto saturation = max(1e-8, curVolVars.saturation(phaseIdx));
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();
210 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
211 partialDerivatives[compIdx][compIdx] += d_storage;
224 template<
class PartialDerivativeMatrix>
226 const Problem& problem,
227 const Element& element,
228 const FVElementGeometry& fvGeometry,
229 const VolumeVariables& curVolVars,
230 const SubControlVolume& scv)
const
235 template<
class PartialDerivativeMatrices,
class T = TypeTag>
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
246 DUNE_THROW(Dune::NotImplemented,
"Analytic flux differentiation only implemented for tpfa");
249 auto rho = [](
const VolumeVariables& volVars)
250 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
253 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
256 static const Scalar upwindWeight = getParam<Scalar>(
"Flux.UpwindWeight");
259 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
260 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
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);
274 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
277 Scalar diffDeriv = 0.0;
280 diffDeriv = useMoles ?
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
286 :
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
289 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
291 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
292 if (!scvf.boundary())
293 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
297 template<
class JacobianMatrix,
class T = TypeTag>
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
309 auto rho = [](
const VolumeVariables& volVars)
310 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
313 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
316 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
319 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
320 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
328 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
330 const auto ti = DiffusionType::calculateTransmissibilities(problem,
335 elemFluxVarsCache[scvf],
337 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
338 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
340 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
342 for (
const auto& scv : scvs(fvGeometry))
345 auto diffDeriv = 0.0;
347 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
348 : ti[compIdx][scv.indexInElement()];
350 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
351 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
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;
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;
365 template<
class PartialDerivativeMatrices>
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
376 curElemVolVars, elemFluxVarsCache, scvf);
379 template<
class PartialDerivativeMatrices>
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
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
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.