3.6-git
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 // checks if the fluid system uses the Richards model index convention
83 static constexpr auto fsCheck = GetPropType<TypeTag, Properties::ModelTraits>::checkFluidSystem(FluidSystem{});
84
85 // phase & component indices
86 static constexpr auto liquidPhaseIdx = FluidSystem::phase0Idx;
87 static constexpr auto gasPhaseIdx = FluidSystem::phase1Idx;
88 static constexpr auto liquidCompIdx = FluidSystem::comp0Idx;
89
90 static constexpr bool enableWaterDiffusionInAir
91 = getPropValue<TypeTag, Properties::EnableWaterDiffusionInAir>();
92
94 struct InvalidElemSol
95 {
96 template<class Index>
97 double operator[] (const Index i) const
98 { static_assert(AlwaysFalse<Index>::value, "Solution-dependent material parameters not supported with analytical differentiation"); return 0.0; }
99 };
100
101public:
102 using ParentType::ParentType;
103
115 NumEqVector computeStorage(const Problem& problem,
116 const SubControlVolume& scv,
117 const VolumeVariables& volVars) const
118 {
119 // partial time derivative of the phase mass
120 NumEqVector storage(0.0);
121 storage[conti0EqIdx] = volVars.porosity()
122 * volVars.density(liquidPhaseIdx)
123 * volVars.saturation(liquidPhaseIdx);
124
125 // for extended Richards we consider water in air
126 if constexpr (enableWaterDiffusionInAir)
127 storage[conti0EqIdx] += volVars.porosity()
128 * volVars.molarDensity(gasPhaseIdx)
129 * volVars.moleFraction(gasPhaseIdx, liquidCompIdx)
130 * FluidSystem::molarMass(liquidCompIdx)
131 * volVars.saturation(gasPhaseIdx);
132
134 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, liquidPhaseIdx);
135 EnergyLocalResidual::fluidPhaseStorage(storage, scv, volVars, gasPhaseIdx);
136 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
137
138 return storage;
139 }
140
141
152 NumEqVector computeFlux(const Problem& problem,
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVolumeVariables& elemVolVars,
156 const SubControlVolumeFace& scvf,
157 const ElementFluxVariablesCache& elemFluxVarsCache) const
158 {
159 FluxVariables fluxVars;
160 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
161
162 NumEqVector flux(0.0);
163 // the physical quantities for which we perform upwinding
164 auto upwindTerm = [](const auto& volVars)
165 { return volVars.density(liquidPhaseIdx)*volVars.mobility(liquidPhaseIdx); };
166
167 flux[conti0EqIdx] = fluxVars.advectiveFlux(liquidPhaseIdx, upwindTerm);
168
169 // for extended Richards we consider water vapor diffusion in air
170 if constexpr (enableWaterDiffusionInAir)
171 {
172 //check for the reference system and adapt units of the diffusive flux accordingly.
173 if (FluxVariables::MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
174 flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx];
175 else
176 flux[conti0EqIdx] += fluxVars.molecularDiffusionFlux(gasPhaseIdx)[liquidCompIdx]*FluidSystem::molarMass(liquidCompIdx);
177 }
178
180 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, liquidPhaseIdx);
181
184 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
185
186 return flux;
187 }
188
199 template<class PartialDerivativeMatrix>
200 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
201 const Problem& problem,
202 const Element& element,
203 const FVElementGeometry& fvGeometry,
204 const VolumeVariables& curVolVars,
205 const SubControlVolume& scv) const
206 {
207 static_assert(!enableWaterDiffusionInAir,
208 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
209 static_assert(!FluidSystem::isCompressible(0),
210 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
211
212 const auto poreVolume = Extrusion::volume(fvGeometry, scv)*curVolVars.porosity()*curVolVars.extrusionFactor();
213 static const auto rho = curVolVars.density(0);
214
215 // partial derivative of storage term w.r.t. p_w
216 // 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
217 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, InvalidElemSol{});
218 partialDerivatives[conti0EqIdx][0] += -rho*poreVolume/this->timeLoop().timeStepSize()*fluidMatrixInteraction.dsw_dpc(curVolVars.capillaryPressure());
219 }
220
233 template<class PartialDerivativeMatrix>
234 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
235 const Problem& problem,
236 const Element& element,
237 const FVElementGeometry& fvGeometry,
238 const VolumeVariables& curVolVars,
239 const SubControlVolume& scv) const
240 { /* TODO maybe forward to problem for the user to implement the source derivatives?*/ }
241
256 template<class PartialDerivativeMatrices, class T = TypeTag>
257 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::cctpfa, void>
258 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
259 const Problem& problem,
260 const Element& element,
261 const FVElementGeometry& fvGeometry,
262 const ElementVolumeVariables& curElemVolVars,
263 const ElementFluxVariablesCache& elemFluxVarsCache,
264 const SubControlVolumeFace& scvf) const
265 {
266 static_assert(!enableWaterDiffusionInAir,
267 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
268 static_assert(!FluidSystem::isCompressible(0),
269 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
270 static_assert(FluidSystem::viscosityIsConstant(0),
271 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
272
273 // get references to the two participating vol vars & parameters
274 const auto insideScvIdx = scvf.insideScvIdx();
275 const auto outsideScvIdx = scvf.outsideScvIdx();
276 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
277 const auto& insideScv = fvGeometry.scv(insideScvIdx);
278 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
279 const auto& insideVolVars = curElemVolVars[insideScvIdx];
280 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
281
282 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
283 static const auto rho = insideVolVars.density(0);
284 static const auto mu = insideVolVars.viscosity(0);
285 static const auto rho_mu = rho/mu;
286
287 // upwind term
288 // evaluate the current wetting phase Darcy flux and resulting upwind weights
290 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
291 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
292 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
293 const auto outsideWeight = 1.0 - insideWeight;
294 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
295
296 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
297 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(outsideElement, outsideScv, InvalidElemSol{});
298
299 // material law derivatives
300 const auto insideSw = insideVolVars.saturation(0);
301 const auto outsideSw = outsideVolVars.saturation(0);
302 const auto insidePc = insideVolVars.capillaryPressure();
303 const auto outsidePc = outsideVolVars.capillaryPressure();
304 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
305 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
306 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
307 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
308
309 // the transmissibility
310 const auto tij = elemFluxVarsCache[scvf].advectionTij();
311
312 // get references to the two participating derivative matrices
313 auto& dI_dI = derivativeMatrices[insideScvIdx];
314 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
315
316 // partial derivative of the wetting phase flux w.r.t. pw
317 dI_dI[conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
318 dI_dJ[conti0EqIdx][0] += -tij*upwindTerm + rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
319 }
320
332 template<class JacobianMatrix, class T = TypeTag>
333 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
334 addFluxDerivatives(JacobianMatrix& A,
335 const Problem& problem,
336 const Element& element,
337 const FVElementGeometry& fvGeometry,
338 const ElementVolumeVariables& curElemVolVars,
339 const ElementFluxVariablesCache& elemFluxVarsCache,
340 const SubControlVolumeFace& scvf) const
341 {
342 static_assert(!enableWaterDiffusionInAir,
343 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
344 static_assert(!FluidSystem::isCompressible(0),
345 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
346 static_assert(FluidSystem::viscosityIsConstant(0),
347 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
348
349 // get references to the two participating vol vars & parameters
350 const auto insideScvIdx = scvf.insideScvIdx();
351 const auto outsideScvIdx = scvf.outsideScvIdx();
352 const auto& insideScv = fvGeometry.scv(insideScvIdx);
353 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
354 const auto& insideVolVars = curElemVolVars[insideScvIdx];
355 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
356
357 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
358 static const auto rho = insideVolVars.density(0);
359 static const auto mu = insideVolVars.viscosity(0);
360 static const auto rho_mu = rho/mu;
361
362 // upwind term
363 // evaluate the current wetting phase Darcy flux and resulting upwind weights
365 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
366 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
367 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
368 const auto outsideWeight = 1.0 - insideWeight;
369 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
370
371 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
372 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, outsideScv, InvalidElemSol{});
373
374 // material law derivatives
375 const auto insideSw = insideVolVars.saturation(0);
376 const auto outsideSw = outsideVolVars.saturation(0);
377 const auto insidePc = insideVolVars.capillaryPressure();
378 const auto outsidePc = outsideVolVars.capillaryPressure();
379 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
380 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
381 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
382 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
383
384 // so far it was the same as for tpfa
385 // the transmissibilities (flux derivatives with respect to all pw-dofs on the element)
386 const auto ti = AdvectionType::calculateTransmissibilities(
387 problem, element, fvGeometry, curElemVolVars, scvf, elemFluxVarsCache[scvf]
388 );
389
390 // get the rows of the jacobian matrix for the inside/outside scv
391 auto& dI_dJ_inside = A[insideScv.dofIndex()];
392 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
393
394 // add the partial derivatives w.r.t all scvs in the element
395 for (const auto& scvJ : scvs(fvGeometry))
396 {
397 // the transmissibily associated with the scvJ
398 const auto& tj = ti[scvJ.indexInElement()];
399 const auto globalJ = scvJ.dofIndex();
400
401 // partial derivative of the wetting phase flux w.r.t. p_w
402 dI_dJ_inside[globalJ][conti0EqIdx][0] += tj*upwindTerm;
403 dI_dJ_outside[globalJ][conti0EqIdx][0] += -tj*upwindTerm;
404 }
405
406 const auto upwindContributionInside = rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
407 const auto upwindContributionOutside = rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
408
409 // additional contribution of the upwind term only for inside and outside dof
410 A[insideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] += upwindContributionInside;
411 A[insideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] += upwindContributionOutside;
412
413 A[outsideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionInside;
414 A[outsideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionOutside;
415 }
416
431 template<class PartialDerivativeMatrices>
432 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
433 const Problem& problem,
434 const Element& element,
435 const FVElementGeometry& fvGeometry,
436 const ElementVolumeVariables& curElemVolVars,
437 const ElementFluxVariablesCache& elemFluxVarsCache,
438 const SubControlVolumeFace& scvf) const
439 {
440 static_assert(!enableWaterDiffusionInAir,
441 "richards/localresidual.hh: Analytic Jacobian not implemented for the water diffusion in air version!");
442 static_assert(!FluidSystem::isCompressible(0),
443 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
444 static_assert(FluidSystem::viscosityIsConstant(0),
445 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
446
447
448 // get references to the two participating vol vars & parameters
449 const auto insideScvIdx = scvf.insideScvIdx();
450 const auto& insideScv = fvGeometry.scv(insideScvIdx);
451 const auto& insideVolVars = curElemVolVars[insideScvIdx];
452 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
453 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
454
455 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
456 static const auto rho = insideVolVars.density(0);
457 static const auto mu = insideVolVars.viscosity(0);
458 static const auto rho_mu = rho/mu;
459
460 // upwind term
461 // evaluate the current wetting phase Darcy flux and resulting upwind weights
463 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
464 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
465 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
466 const auto outsideWeight = 1.0 - insideWeight;
467 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
468
469 // material law derivatives
470 const auto insideSw = insideVolVars.saturation(0);
471 const auto insidePc = insideVolVars.capillaryPressure();
472 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
473 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
474
475 // the transmissibility
476 const auto tij = elemFluxVarsCache[scvf].advectionTij();
477
478 // partial derivative of the wetting phase flux w.r.t. pw
479 derivativeMatrices[insideScvIdx][conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
480 }
481
493 template<class PartialDerivativeMatrices>
494 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
495 const Problem& problem,
496 const Element& element,
497 const FVElementGeometry& fvGeometry,
498 const ElementVolumeVariables& curElemVolVars,
499 const ElementFluxVariablesCache& elemFluxVarsCache,
500 const SubControlVolumeFace& scvf) const
501 {
502 if constexpr(Detail::hasAddRobinFluxDerivatives<Problem,
503 PartialDerivativeMatrices&, Element, FVElementGeometry,
504 ElementVolumeVariables, ElementFluxVariablesCache, SubControlVolumeFace>()
505 )
506 problem.addRobinFluxDerivatives(derivativeMatrices, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
507 }
508
509private:
510 Implementation *asImp_()
511 { return static_cast<Implementation *> (this); }
512
513 const Implementation *asImp_() const
514 { return static_cast<const Implementation *> (this); }
515};
516
517} // end namespace Dumux
518
519#endif
Type traits.
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Distance implementation details.
Definition: cvfelocalresidual.hh:37
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
constexpr CCTpfa cctpfa
Definition: method.hh:134
constexpr Box box
Definition: method.hh:136
Template which always yields a false value.
Definition: typetraits.hh:35
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:234
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::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:258
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:432
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:115
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:200
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:152
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:494
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::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:334
Declares all properties used in Dumux.