26#ifndef DUMUX_TRACER_LOCAL_RESIDUAL_HH
27#define DUMUX_TRACER_LOCAL_RESIDUAL_HH
41template<
class TypeTag>
48 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
49 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
54 using Element =
typename GridView::template Codim<0>::Entity;
60 static constexpr bool useMoles = getPropValue<TypeTag, Properties::UseMoles>();
61 static constexpr int phaseIdx = 0;
64 using ParentType::ParentType;
78 const SubControlVolume& scv,
79 const VolumeVariables& volVars)
const
81 NumEqVector storage(0.0);
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);
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);
118 const Element& element,
119 const FVElementGeometry& fvGeometry,
120 const ElementVolumeVariables& elemVolVars,
121 const SubControlVolumeFace& scvf,
122 const ElementFluxVariablesCache& elemFluxVarsCache)
const
124 FluxVariables fluxVars;
125 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
127 NumEqVector flux(0.0);
128 const auto diffusiveFluxes = fluxVars.molecularDiffusionFlux(phaseIdx);
130 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
134 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
137 auto upwindTerm = [compIdx](
const auto& volVars)
138 {
return volVars.molarDensity()*volVars.moleFraction(phaseIdx, compIdx); };
141 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
144 flux[compIdx] += diffusiveFluxes[compIdx]/FluidSystem::molarMass(compIdx);
146 flux[compIdx] += diffusiveFluxes[compIdx];
148 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
154 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
157 auto upwindTerm = [compIdx](
const auto& volVars)
158 {
return volVars.density()*volVars.massFraction(phaseIdx, compIdx); };
161 flux[compIdx] += fluxVars.advectiveFlux(phaseIdx, upwindTerm);
164 flux[compIdx] += diffusiveFluxes[compIdx];
166 flux[compIdx] += diffusiveFluxes[compIdx]*FluidSystem::molarMass(compIdx);
168 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
185 template<
class PartialDerivativeMatrix>
187 const Problem& problem,
188 const Element& element,
189 const FVElementGeometry& fvGeometry,
190 const VolumeVariables& curVolVars,
191 const SubControlVolume& scv)
const
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();
197 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
198 partialDerivatives[compIdx][compIdx] += d_storage;
211 template<
class PartialDerivativeMatrix>
213 const Problem& problem,
214 const Element& element,
215 const FVElementGeometry& fvGeometry,
216 const VolumeVariables& curVolVars,
217 const SubControlVolume& scv)
const
222 template<
class PartialDerivativeMatrices,
class T = TypeTag>
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
233 auto rho = [](
const VolumeVariables& volVars)
234 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
237 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
240 static const Scalar upwindWeight = getParam<Scalar>(
"Flux.UpwindWeight");
243 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
244 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
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);
258 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
261 Scalar diffDeriv = 0.0;
264 diffDeriv = useMoles ?
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)/FluidSystem::molarMass(compIdx)
270 :
massOrMolarDensity*fluxCache.diffusionTij(phaseIdx, compIdx)*FluidSystem::molarMass(compIdx);
273 DUNE_THROW(Dune::NotImplemented,
"other reference systems than mass and molar averaged are not implemented");
275 derivativeMatrices[scvf.insideScvIdx()][compIdx][compIdx] += (advDerivII + diffDeriv);
276 derivativeMatrices[scvf.outsideScvIdx()][compIdx][compIdx] += (advDerivIJ - diffDeriv);
280 template<
class JacobianMatrix,
class T = TypeTag>
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
292 auto rho = [](
const VolumeVariables& volVars)
293 {
return useMoles ? volVars.molarDensity() : volVars.density(); };
296 const auto volFlux = problem.spatialParams().volumeFlux(element, fvGeometry, curElemVolVars, scvf);
299 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
302 const auto& insideVolVars = curElemVolVars[scvf.insideScvIdx()];
303 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
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;
311 static constexpr auto referenceSystemFormulation = FluxVariables::MolecularDiffusionType::referenceSystemFormulation();
313 const auto ti = DiffusionType::calculateTransmissibilities(problem,
318 elemFluxVarsCache[scvf],
320 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
321 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
323 for (
int compIdx = 0; compIdx < numComponents; ++compIdx)
325 for (
const auto& scv : scvs(fvGeometry))
328 auto diffDeriv = 0.0;
330 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]/FluidSystem::molarMass(compIdx)
331 : ti[compIdx][scv.indexInElement()];
333 diffDeriv += useMoles ? ti[compIdx][scv.indexInElement()]
334 : ti[compIdx][scv.indexInElement()]*FluidSystem::molarMass(compIdx);
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;
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;
348 template<
class PartialDerivativeMatrices>
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
359 curElemVolVars, elemFluxVarsCache, scvf);
362 template<
class PartialDerivativeMatrices>
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
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.