3.3.0
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
35
37
38namespace Dumux {
39
45template<class TypeTag>
46class RichardsLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
47{
49
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
57 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
59 using FVElementGeometry = typename GridGeometry::LocalView;
60 using Extrusion = Extrusion_t<GridGeometry>;
61 using SubControlVolume = typename GridGeometry::SubControlVolume;
62 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
63 using GridView = typename GridGeometry::GridView;
64 using Element = typename GridView::template Codim<0>::Entity;
68 // first index for the mass balance
69 enum { conti0EqIdx = Indices::conti0EqIdx };
70
71 // phase indices
72 enum {
73 liquidPhaseIdx = FluidSystem::liquidPhaseIdx,
74 gasPhaseIdx = FluidSystem::gasPhaseIdx,
75 liquidCompIdx = FluidSystem::liquidCompIdx
76 };
77
78 static constexpr bool enableWaterDiffusionInAir
79 = getPropValue<TypeTag, Properties::EnableWaterDiffusionInAir>();
80
82 struct InvalidElemSol
83 {
84 template<class Index>
85 double operator[] (const Index i) const
86 { static_assert(AlwaysFalse<Index>::value, "Solution-dependent material parameters not supported with analytical differentiation"); return 0.0; }
87 };
88
89public:
90 using ParentType::ParentType;
91
103 NumEqVector computeStorage(const Problem& problem,
104 const SubControlVolume& scv,
105 const VolumeVariables& volVars) const
106 {
107 // partial time derivative of the phase mass
108 NumEqVector storage(0.0);
109 storage[conti0EqIdx] = volVars.porosity()
110 * volVars.density(liquidPhaseIdx)
111 * volVars.saturation(liquidPhaseIdx);
112
113 // for extended Richards we consider water in air
114 if (enableWaterDiffusionInAir)
115 storage[conti0EqIdx] += volVars.porosity()
116 * volVars.molarDensity(gasPhaseIdx)
117 * volVars.moleFraction(gasPhaseIdx, liquidCompIdx)
118 * FluidSystem::molarMass(liquidCompIdx)
119 * volVars.saturation(gasPhaseIdx);
120
122 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, liquidPhaseIdx);
123 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, gasPhaseIdx);
124 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
125
126 return storage;
127 }
128
129
140 NumEqVector computeFlux(const Problem& problem,
141 const Element& element,
142 const FVElementGeometry& fvGeometry,
143 const ElementVolumeVariables& elemVolVars,
144 const SubControlVolumeFace& scvf,
145 const ElementFluxVariablesCache& elemFluxVarsCache) const
146 {
147 FluxVariables fluxVars;
148 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
149
150 NumEqVector flux(0.0);
151 // the physical quantities for which we perform upwinding
152 auto upwindTerm = [](const auto& volVars)
153 { return volVars.density(liquidPhaseIdx)*volVars.mobility(liquidPhaseIdx); };
154
155 flux[conti0EqIdx] = fluxVars.advectiveFlux(liquidPhaseIdx, upwindTerm);
156
157 // for extended Richards we consider water vapor diffusion in air
158 if (enableWaterDiffusionInAir)
159 {
160 //check for the reference system and adapt units of the diffusive flux accordingly.
161 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
162 flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx];
163 else
164 flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx]*FluidSystem::molarMass(liquidCompIdx);
165 }
166
168 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, liquidPhaseIdx);
169
172 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
173
174 return flux;
175 }
176
187 template<class PartialDerivativeMatrix>
188 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
189 const Problem& problem,
190 const Element& element,
191 const FVElementGeometry& fvGeometry,
192 const VolumeVariables& curVolVars,
193 const SubControlVolume& scv) const
194 {
195 static_assert(!enableWaterDiffusionInAir,
196 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
197 static_assert(!FluidSystem::isCompressible(0),
198 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
199
200 const auto poreVolume = Extrusion::volume(scv)*curVolVars.porosity();
201 static const auto rho = curVolVars.density(0);
202
203 // old material law interface is deprecated: Replace this by
204 // const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
205 // after the release of 3.3, when the deprecated interface is no longer supported
206 const auto fluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element, scv, InvalidElemSol{});
207
208 // partial derivative of storage term w.r.t. p_w
209 // 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
210 partialDerivatives[conti0EqIdx][0] += -rho*poreVolume/this->timeLoop().timeStepSize()*fluidMatrixInteraction.dsw_dpc(curVolVars.capillaryPressure());
211 }
212
225 template<class PartialDerivativeMatrix>
226 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
227 const Problem& problem,
228 const Element& element,
229 const FVElementGeometry& fvGeometry,
230 const VolumeVariables& curVolVars,
231 const SubControlVolume& scv) const
232 { /* TODO maybe forward to problem for the user to implement the source derivatives?*/ }
233
248 template<class PartialDerivativeMatrices, class T = TypeTag>
249 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::cctpfa, void>
250 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
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
257 {
258 static_assert(!enableWaterDiffusionInAir,
259 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
260 static_assert(!FluidSystem::isCompressible(0),
261 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
262 static_assert(!FluidSystem::viscosityIsConstant(0),
263 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
264
265 // get references to the two participating vol vars & parameters
266 const auto insideScvIdx = scvf.insideScvIdx();
267 const auto outsideScvIdx = scvf.outsideScvIdx();
268 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
269 const auto& insideScv = fvGeometry.scv(insideScvIdx);
270 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
271 const auto& insideVolVars = curElemVolVars[insideScvIdx];
272 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
273
274 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
275 static const auto rho = insideVolVars.density(0);
276 static const auto mu = insideVolVars.viscosity(0);
277 static const auto rho_mu = rho/mu;
278
279 // upwind term
280 // evaluate the current wetting phase Darcy flux and resulting upwind weights
282 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
283 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
284 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
285 const auto outsideWeight = 1.0 - insideWeight;
286 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
287
288 // old material law interface is deprecated: Replace this by
289 // const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
290 // after the release of 3.3, when the deprecated interface is no longer supported
291 const auto insideFluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element, insideScv, InvalidElemSol{});
292 const auto outsideFluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), outsideElement, outsideScv, InvalidElemSol{});
293
294 // material law derivatives
295 const auto insideSw = insideVolVars.saturation(0);
296 const auto outsideSw = outsideVolVars.saturation(0);
297 const auto insidePc = insideVolVars.capillaryPressure();
298 const auto outsidePc = outsideVolVars.capillaryPressure();
299 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
300 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
301 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
302 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
303
304 // the transmissibility
305 const auto tij = elemFluxVarsCache[scvf].advectionTij();
306
307 // get references to the two participating derivative matrices
308 auto& dI_dI = derivativeMatrices[insideScvIdx];
309 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
310
311 // partial derivative of the wetting phase flux w.r.t. pw
312 dI_dI[conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
313 dI_dJ[conti0EqIdx][0] += -tij*upwindTerm + rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
314 }
315
327 template<class JacobianMatrix, class T = TypeTag>
328 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethod::box, void>
329 addFluxDerivatives(JacobianMatrix& A,
330 const Problem& problem,
331 const Element& element,
332 const FVElementGeometry& fvGeometry,
333 const ElementVolumeVariables& curElemVolVars,
334 const ElementFluxVariablesCache& elemFluxVarsCache,
335 const SubControlVolumeFace& scvf) const
336 {
337 static_assert(!enableWaterDiffusionInAir,
338 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
339 static_assert(!FluidSystem::isCompressible(0),
340 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
341 static_assert(!FluidSystem::viscosityIsConstant(0),
342 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
343
344 // get references to the two participating vol vars & parameters
345 const auto insideScvIdx = scvf.insideScvIdx();
346 const auto outsideScvIdx = scvf.outsideScvIdx();
347 const auto& insideScv = fvGeometry.scv(insideScvIdx);
348 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
349 const auto& insideVolVars = curElemVolVars[insideScvIdx];
350 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
351
352 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
353 static const auto rho = insideVolVars.density(0);
354 static const auto mu = insideVolVars.viscosity(0);
355 static const auto rho_mu = rho/mu;
356
357 // upwind term
358 // evaluate the current wetting phase Darcy flux and resulting upwind weights
360 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
361 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
362 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
363 const auto outsideWeight = 1.0 - insideWeight;
364 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
365
366 // old material law interface is deprecated: Replace this by
367 // const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
368 // after the release of 3.3, when the deprecated interface is no longer supported
369 const auto insideFluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element, insideScv, InvalidElemSol{});
370 const auto outsideFluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element, outsideScv, InvalidElemSol{});
371
372 // material law derivatives
373 const auto insideSw = insideVolVars.saturation(0);
374 const auto outsideSw = outsideVolVars.saturation(0);
375 const auto insidePc = insideVolVars.capillaryPressure();
376 const auto outsidePc = outsideVolVars.capillaryPressure();
377 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
378 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
379 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
380 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
381
382 // so far it was the same as for tpfa
383 // the transmissibilities (flux derivatives with respect to all pw-dofs on the element)
384 const auto ti = AdvectionType::calculateTransmissibilities(
385 problem, element, fvGeometry, curElemVolVars, scvf, elemFluxVarsCache[scvf]
386 );
387
388 // get the rows of the jacobian matrix for the inside/outside scv
389 auto& dI_dJ_inside = A[insideScv.dofIndex()];
390 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
391
392 // add the partial derivatives w.r.t all scvs in the element
393 for (const auto& scvJ : scvs(fvGeometry))
394 {
395 // the transmissibily associated with the scvJ
396 const auto& tj = ti[scvJ.indexInElement()];
397 const auto globalJ = scvJ.dofIndex();
398
399 // partial derivative of the wetting phase flux w.r.t. p_w
400 dI_dJ_inside[globalJ][conti0EqIdx][0] += tj*upwindTerm;
401 dI_dJ_outside[globalJ][conti0EqIdx][0] += -tj*upwindTerm;
402 }
403
404 const auto upwindContributionInside = rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
405 const auto upwindContributionOutside = rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
406
407 // additional constribution of the upwind term only for inside and outside dof
408 A[insideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] += upwindContributionInside;
409 A[insideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] += upwindContributionOutside;
410
411 A[outsideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionInside;
412 A[outsideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionOutside;
413 }
414
429 template<class PartialDerivativeMatrices>
430 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
431 const Problem& problem,
432 const Element& element,
433 const FVElementGeometry& fvGeometry,
434 const ElementVolumeVariables& curElemVolVars,
435 const ElementFluxVariablesCache& elemFluxVarsCache,
436 const SubControlVolumeFace& scvf) const
437 {
438 static_assert(!enableWaterDiffusionInAir,
439 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
440 static_assert(!FluidSystem::isCompressible(0),
441 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
442 static_assert(!FluidSystem::viscosityIsConstant(0),
443 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
444
445
446 // get references to the two participating vol vars & parameters
447 const auto insideScvIdx = scvf.insideScvIdx();
448 const auto& insideScv = fvGeometry.scv(insideScvIdx);
449 const auto& insideVolVars = curElemVolVars[insideScvIdx];
450 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
451
452 // old material law interface is deprecated: Replace this by
453 // const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, elemSol);
454 // after the release of 3.3, when the deprecated interface is no longer supported
455 const auto insideFluidMatrixInteraction = Deprecated::makePcKrSw(Scalar{}, problem.spatialParams(), element, insideScv, InvalidElemSol{});
456
457 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
458 static const auto rho = insideVolVars.density(0);
459 static const auto mu = insideVolVars.viscosity(0);
460 static const auto rho_mu = rho/mu;
461
462 // upwind term
463 // evaluate the current wetting phase Darcy flux and resulting upwind weights
465 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
466 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
467 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
468 const auto outsideWeight = 1.0 - insideWeight;
469 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
470
471 // material law derivatives
472 const auto insideSw = insideVolVars.saturation(0);
473 const auto insidePc = insideVolVars.capillaryPressure();
474 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
475 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
476
477 // the transmissibility
478 const auto tij = elemFluxVarsCache[scvf].advectionTij();
479
480 // partial derivative of the wetting phase flux w.r.t. pw
481 derivativeMatrices[insideScvIdx][conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
482 }
483
495 template<class PartialDerivativeMatrices>
496 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
497 const Problem& problem,
498 const Element& element,
499 const FVElementGeometry& fvGeometry,
500 const ElementVolumeVariables& curElemVolVars,
501 const ElementFluxVariablesCache& elemFluxVarsCache,
502 const SubControlVolumeFace& scvf) const
503 { /* TODO maybe forward to problem for the user to implement the Robin derivatives?*/ }
504
505private:
506 Implementation *asImp_()
507 { return static_cast<Implementation *> (this); }
508
509 const Implementation *asImp_() const
510 { return static_cast<const Implementation *> (this); }
511};
512
513} // end namespace Dumux
514
515#endif
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Helpers for deprecation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
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 (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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:47
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:226
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:430
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:103
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:188
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:140
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:250
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:496
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:329
Declares all properties used in Dumux.