3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
porousmediumflow/richards/localresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
26#ifndef DUMUX_RICHARDS_LOCAL_RESIDUAL_HH
27#define DUMUX_RICHARDS_LOCAL_RESIDUAL_HH
28
36
37namespace Dumux::Detail {
38// helper structs and functions detecting if the user-defined problem class
39// implements addRobinFluxDerivatives
40template <typename T, typename ...Ts>
41using RobinDerivDetector = decltype(std::declval<T>().addRobinFluxDerivatives(std::declval<Ts>()...));
42
43template<class T, typename ...Args>
44static constexpr bool hasAddRobinFluxDerivatives()
45{ return Dune::Std::is_detected<RobinDerivDetector, T, Args...>::value; }
46
47} // end namespace Dumux::Detail
48
49namespace Dumux {
50
56template<class TypeTag>
57class RichardsLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
58{
60
66 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
68 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
70 using FVElementGeometry = typename GridGeometry::LocalView;
71 using Extrusion = Extrusion_t<GridGeometry>;
72 using SubControlVolume = typename GridGeometry::SubControlVolume;
73 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
74 using GridView = typename GridGeometry::GridView;
75 using Element = typename GridView::template Codim<0>::Entity;
79 // first index for the mass balance
80 enum { conti0EqIdx = Indices::conti0EqIdx };
81
82 // phase indices
83 enum {
84 liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
85 gasPhaseIdx = FluidSystem::gasPhaseIdx,
86 liquidCompIdx = FluidSystem::liquidCompIdx
87 };
88
89 static constexpr bool enableWaterDiffusionInAir
90 = getPropValue<TypeTag, Properties::EnableWaterDiffusionInAir>();
91
93 struct InvalidElemSol
94 {
95 template<class Index>
96 double operator[] (const Index i) const
97 { static_assert(AlwaysFalse<Index>::value, "Solution-dependent material parameters not supported with analytical differentiation"); return 0.0; }
98 };
99
100public:
101 using ParentType::ParentType;
102
114 NumEqVector computeStorage(const Problem& problem,
115 const SubControlVolume& scv,
116 const VolumeVariables& volVars) const
117 {
118 // partial time derivative of the phase mass
119 NumEqVector storage(0.0);
120 storage[conti0EqIdx] = volVars.porosity()
121 * volVars.density(liquidPhaseIdx)
122 * volVars.saturation(liquidPhaseIdx);
123
124 // for extended Richards we consider water in air
125 if constexpr (enableWaterDiffusionInAir)
126 storage[conti0EqIdx] += volVars.porosity()
127 * volVars.molarDensity(gasPhaseIdx)
128 * volVars.moleFraction(gasPhaseIdx, liquidCompIdx)
129 * FluidSystem::molarMass(liquidCompIdx)
130 * volVars.saturation(gasPhaseIdx);
131
133 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, liquidPhaseIdx);
134 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, gasPhaseIdx);
135 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
136
137 return storage;
138 }
139
140
151 NumEqVector computeFlux(const Problem& problem,
152 const Element& element,
153 const FVElementGeometry& fvGeometry,
154 const ElementVolumeVariables& elemVolVars,
155 const SubControlVolumeFace& scvf,
156 const ElementFluxVariablesCache& elemFluxVarsCache) const
157 {
158 FluxVariables fluxVars;
159 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
160
161 NumEqVector flux(0.0);
162 // the physical quantities for which we perform upwinding
163 auto upwindTerm = [](const auto& volVars)
164 { return volVars.density(liquidPhaseIdx)*volVars.mobility(liquidPhaseIdx); };
165
166 flux[conti0EqIdx] = fluxVars.advectiveFlux(liquidPhaseIdx, upwindTerm);
167
168 // for extended Richards we consider water vapor diffusion in air
169 if constexpr (enableWaterDiffusionInAir)
170 {
171 //check for the reference system and adapt units of the diffusive flux accordingly.
172 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
173 flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx];
174 else
175 flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx]*FluidSystem::molarMass(liquidCompIdx);
176 }
177
179 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, liquidPhaseIdx);
180
183 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
184
185 return flux;
186 }
187
198 template<class PartialDerivativeMatrix>
199 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
200 const Problem& problem,
201 const Element& element,
202 const FVElementGeometry& fvGeometry,
203 const VolumeVariables& curVolVars,
204 const SubControlVolume& scv) const
205 {
206 static_assert(!enableWaterDiffusionInAir,
207 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
208 static_assert(!FluidSystem::isCompressible(0),
209 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
210
211 const auto poreVolume = Extrusion::volume(scv)*curVolVars.porosity()*curVolVars.extrusionFactor();
212 static const auto rho = curVolVars.density(0);
213
214 // partial derivative of storage term w.r.t. p_w
215 // d(Sw*rho*phi*V/dt)/dpw = rho*phi*V/dt*dsw/dpw = rho*phi*V/dt*dsw/dpc*dpc/dpw = -rho*phi*V/dt*dsw/dpc
216 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, InvalidElemSol{});
217 partialDerivatives[conti0EqIdx][0] += -rho*poreVolume/this->timeLoop().timeStepSize()*fluidMatrixInteraction.dsw_dpc(curVolVars.capillaryPressure());
218 }
219
232 template<class PartialDerivativeMatrix>
233 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
234 const Problem& problem,
235 const Element& element,
236 const FVElementGeometry& fvGeometry,
237 const VolumeVariables& curVolVars,
238 const SubControlVolume& scv) const
239 { /* TODO maybe forward to problem for the user to implement the source derivatives?*/ }
240
255 template<class PartialDerivativeMatrices, class T = TypeTag>
256 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::cctpfa, void>
257 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
258 const Problem& problem,
259 const Element& element,
260 const FVElementGeometry& fvGeometry,
261 const ElementVolumeVariables& curElemVolVars,
262 const ElementFluxVariablesCache& elemFluxVarsCache,
263 const SubControlVolumeFace& scvf) const
264 {
265 static_assert(!enableWaterDiffusionInAir,
266 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
267 static_assert(!FluidSystem::isCompressible(0),
268 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
269 static_assert(!FluidSystem::viscosityIsConstant(0),
270 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
271
272 // get references to the two participating vol vars & parameters
273 const auto insideScvIdx = scvf.insideScvIdx();
274 const auto outsideScvIdx = scvf.outsideScvIdx();
275 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
276 const auto& insideScv = fvGeometry.scv(insideScvIdx);
277 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
278 const auto& insideVolVars = curElemVolVars[insideScvIdx];
279 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
280
281 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
282 static const auto rho = insideVolVars.density(0);
283 static const auto mu = insideVolVars.viscosity(0);
284 static const auto rho_mu = rho/mu;
285
286 // upwind term
287 // evaluate the current wetting phase Darcy flux and resulting upwind weights
289 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
290 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
291 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
292 const auto outsideWeight = 1.0 - insideWeight;
293 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
294
295 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
296 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(outsideElement, outsideScv, InvalidElemSol{});
297
298 // material law derivatives
299 const auto insideSw = insideVolVars.saturation(0);
300 const auto outsideSw = outsideVolVars.saturation(0);
301 const auto insidePc = insideVolVars.capillaryPressure();
302 const auto outsidePc = outsideVolVars.capillaryPressure();
303 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
304 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
305 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
306 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
307
308 // the transmissibility
309 const auto tij = elemFluxVarsCache[scvf].advectionTij();
310
311 // get references to the two participating derivative matrices
312 auto& dI_dI = derivativeMatrices[insideScvIdx];
313 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
314
315 // partial derivative of the wetting phase flux w.r.t. pw
316 dI_dI[conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
317 dI_dJ[conti0EqIdx][0] += -tij*upwindTerm + rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
318 }
319
331 template<class JacobianMatrix, class T = TypeTag>
332 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::box, void>
333 addFluxDerivatives(JacobianMatrix& A,
334 const Problem& problem,
335 const Element& element,
336 const FVElementGeometry& fvGeometry,
337 const ElementVolumeVariables& curElemVolVars,
338 const ElementFluxVariablesCache& elemFluxVarsCache,
339 const SubControlVolumeFace& scvf) const
340 {
341 static_assert(!enableWaterDiffusionInAir,
342 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
343 static_assert(!FluidSystem::isCompressible(0),
344 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
345 static_assert(!FluidSystem::viscosityIsConstant(0),
346 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
347
348 // get references to the two participating vol vars & parameters
349 const auto insideScvIdx = scvf.insideScvIdx();
350 const auto outsideScvIdx = scvf.outsideScvIdx();
351 const auto& insideScv = fvGeometry.scv(insideScvIdx);
352 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
353 const auto& insideVolVars = curElemVolVars[insideScvIdx];
354 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
355
356 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
357 static const auto rho = insideVolVars.density(0);
358 static const auto mu = insideVolVars.viscosity(0);
359 static const auto rho_mu = rho/mu;
360
361 // upwind term
362 // evaluate the current wetting phase Darcy flux and resulting upwind weights
364 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
365 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
366 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
367 const auto outsideWeight = 1.0 - insideWeight;
368 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
369
370 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
371 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, outsideScv, InvalidElemSol{});
372
373 // material law derivatives
374 const auto insideSw = insideVolVars.saturation(0);
375 const auto outsideSw = outsideVolVars.saturation(0);
376 const auto insidePc = insideVolVars.capillaryPressure();
377 const auto outsidePc = outsideVolVars.capillaryPressure();
378 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
379 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
380 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
381 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
382
383 // so far it was the same as for tpfa
384 // the transmissibilities (flux derivatives with respect to all pw-dofs on the element)
385 const auto ti = AdvectionType::calculateTransmissibilities(
386 problem, element, fvGeometry, curElemVolVars, scvf, elemFluxVarsCache[scvf]
387 );
388
389 // get the rows of the jacobian matrix for the inside/outside scv
390 auto& dI_dJ_inside = A[insideScv.dofIndex()];
391 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
392
393 // add the partial derivatives w.r.t all scvs in the element
394 for (const auto& scvJ : scvs(fvGeometry))
395 {
396 // the transmissibily associated with the scvJ
397 const auto& tj = ti[scvJ.indexInElement()];
398 const auto globalJ = scvJ.dofIndex();
399
400 // partial derivative of the wetting phase flux w.r.t. p_w
401 dI_dJ_inside[globalJ][conti0EqIdx][0] += tj*upwindTerm;
402 dI_dJ_outside[globalJ][conti0EqIdx][0] += -tj*upwindTerm;
403 }
404
405 const auto upwindContributionInside = rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
406 const auto upwindContributionOutside = rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
407
408 // additional constribution of the upwind term only for inside and outside dof
409 A[insideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] += upwindContributionInside;
410 A[insideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] += upwindContributionOutside;
411
412 A[outsideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionInside;
413 A[outsideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionOutside;
414 }
415
430 template<class PartialDerivativeMatrices>
431 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
432 const Problem& problem,
433 const Element& element,
434 const FVElementGeometry& fvGeometry,
435 const ElementVolumeVariables& curElemVolVars,
436 const ElementFluxVariablesCache& elemFluxVarsCache,
437 const SubControlVolumeFace& scvf) const
438 {
439 static_assert(!enableWaterDiffusionInAir,
440 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
441 static_assert(!FluidSystem::isCompressible(0),
442 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
443 static_assert(!FluidSystem::viscosityIsConstant(0),
444 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
445
446
447 // get references to the two participating vol vars & parameters
448 const auto insideScvIdx = scvf.insideScvIdx();
449 const auto& insideScv = fvGeometry.scv(insideScvIdx);
450 const auto& insideVolVars = curElemVolVars[insideScvIdx];
451 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
452 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
453
454 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
455 static const auto rho = insideVolVars.density(0);
456 static const auto mu = insideVolVars.viscosity(0);
457 static const auto rho_mu = rho/mu;
458
459 // upwind term
460 // evaluate the current wetting phase Darcy flux and resulting upwind weights
462 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
463 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
464 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
465 const auto outsideWeight = 1.0 - insideWeight;
466 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
467
468 // material law derivatives
469 const auto insideSw = insideVolVars.saturation(0);
470 const auto insidePc = insideVolVars.capillaryPressure();
471 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
472 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
473
474 // the transmissibility
475 const auto tij = elemFluxVarsCache[scvf].advectionTij();
476
477 // partial derivative of the wetting phase flux w.r.t. pw
478 derivativeMatrices[insideScvIdx][conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
479 }
480
492 template<class PartialDerivativeMatrices>
493 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
494 const Problem& problem,
495 const Element& element,
496 const FVElementGeometry& fvGeometry,
497 const ElementVolumeVariables& curElemVolVars,
498 const ElementFluxVariablesCache& elemFluxVarsCache,
499 const SubControlVolumeFace& scvf) const
500 {
501 if constexpr(Detail::hasAddRobinFluxDerivatives<Problem,
502 PartialDerivativeMatrices&, Element, FVElementGeometry,
503 ElementVolumeVariables, ElementFluxVariablesCache, SubControlVolumeFace>()
504 )
505 problem.addRobinFluxDerivatives(derivativeMatrices, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
506 }
507
508private:
509 Implementation *asImp_()
510 { return static_cast<Implementation *> (this); }
511
512 const Implementation *asImp_() const
513 { return static_cast<const Implementation *> (this); }
514};
515
516} // end namespace Dumux
517
518#endif
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...
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
Definition: adapt.hh:29
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
Definition: gridcapabilities.hh:60
decltype(std::declval< T >().addRobinFluxDerivatives(std::declval< Ts >()...)) RobinDerivDetector
Definition: porousmediumflow/richards/localresidual.hh:41
static constexpr bool hasAddRobinFluxDerivatives()
Definition: porousmediumflow/richards/localresidual.hh:44
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Template which always yields a false value.
Definition: typetraits.hh:36
Element-wise calculation of the Jacobian matrix for problems using the Richards fully implicit models...
Definition: porousmediumflow/richards/localresidual.hh:58
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Adds source derivatives for wetting and nonwetting phase.
Definition: porousmediumflow/richards/localresidual.hh:233
void addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Adds cell-centered Dirichlet flux derivatives for wetting and nonwetting phase.
Definition: porousmediumflow/richards/localresidual.hh:431
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluates the rate of change of all conservation quantites (e.g. phase mass) within a sub-control vol...
Definition: porousmediumflow/richards/localresidual.hh:114
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Adds the storage derivative.
Definition: porousmediumflow/richards/localresidual.hh:199
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluates the mass flux over a face of a sub control volume.
Definition: porousmediumflow/richards/localresidual.hh:151
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::cctpfa, void > addFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Adds flux derivatives for wetting and nonwetting phase for cell-centered FVM using TPFA.
Definition: porousmediumflow/richards/localresidual.hh:257
void addRobinFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Adds Robin flux derivatives for wetting and nonwetting phase.
Definition: porousmediumflow/richards/localresidual.hh:493
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
Adds flux derivatives for box method.
Definition: porousmediumflow/richards/localresidual.hh:333
Declares all properties used in Dumux.