version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_RICHARDS_LOCAL_RESIDUAL_HH
15#define DUMUX_RICHARDS_LOCAL_RESIDUAL_HH
16
24
25namespace Dumux::Detail {
26// helper structs and functions detecting if the user-defined problem class
27// implements addRobinFluxDerivatives
28template <typename T, typename ...Ts>
29using RobinDerivDetector = decltype(std::declval<T>().addRobinFluxDerivatives(std::declval<Ts>()...));
30
31template<class T, typename ...Args>
32static constexpr bool hasAddRobinFluxDerivatives()
33{ return Dune::Std::is_detected<RobinDerivDetector, T, Args...>::value; }
34
35} // end namespace Dumux::Detail
36
37namespace Dumux {
38
44template<class TypeTag>
45class RichardsLocalResidual : public GetPropType<TypeTag, Properties::BaseLocalResidual>
46{
48
54 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
56 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
58 using FVElementGeometry = typename GridGeometry::LocalView;
59 using Extrusion = Extrusion_t<GridGeometry>;
60 using SubControlVolume = typename GridGeometry::SubControlVolume;
61 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
62 using GridView = typename GridGeometry::GridView;
63 using Element = typename GridView::template Codim<0>::Entity;
67 // first index for the mass balance
68 enum { conti0EqIdx = Indices::conti0EqIdx };
69
70 // checks if the fluid system uses the Richards model index convention
71 static constexpr auto fsCheck = GetPropType<TypeTag, Properties::ModelTraits>::checkFluidSystem(FluidSystem{});
72
73 // phase & component indices
74 static constexpr auto liquidPhaseIdx = FluidSystem::phase0Idx;
75 static constexpr auto gasPhaseIdx = FluidSystem::phase1Idx;
76 static constexpr auto liquidCompIdx = FluidSystem::comp0Idx;
77
79 struct InvalidElemSol
80 {
81 template<class Index>
82 double operator[] (const Index i) const
83 { static_assert(AlwaysFalse<Index>::value, "Solution-dependent material parameters not supported with analytical differentiation"); return 0.0; }
84 };
85
86public:
87 using ParentType::ParentType;
88
100 NumEqVector computeStorage(const Problem& problem,
101 const SubControlVolume& scv,
102 const VolumeVariables& volVars) const
103 {
104 // partial time derivative of the phase mass
105 NumEqVector storage(0.0);
106 storage[conti0EqIdx] = volVars.porosity()
107 * volVars.density(liquidPhaseIdx)
108 * volVars.saturation(liquidPhaseIdx);
109
111 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, liquidPhaseIdx);
112 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, gasPhaseIdx);
113 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
114
115 return storage;
116 }
117
118
129 NumEqVector computeFlux(const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache) const
135 {
136 FluxVariables fluxVars;
137 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
138
139 NumEqVector flux(0.0);
140 // the physical quantities for which we perform upwinding
141 auto upwindTerm = [](const auto& volVars)
142 { return volVars.density(liquidPhaseIdx)*volVars.mobility(liquidPhaseIdx); };
143
144 flux[conti0EqIdx] = fluxVars.advectiveFlux(liquidPhaseIdx, upwindTerm);
145
147 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, liquidPhaseIdx);
148
151 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
152
153 return flux;
154 }
155
166 template<class PartialDerivativeMatrix>
167 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
168 const Problem& problem,
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const VolumeVariables& curVolVars,
172 const SubControlVolume& scv) const
173 {
174 static_assert(!FluidSystem::isCompressible(0),
175 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
176
177 const auto poreVolume = Extrusion::volume(fvGeometry, scv)*curVolVars.porosity()*curVolVars.extrusionFactor();
178 static const auto rho = curVolVars.density(0);
179
180 // partial derivative of storage term w.r.t. p_w
181 // 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
182 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, InvalidElemSol{});
183 partialDerivatives[conti0EqIdx][0] += -rho*poreVolume/this->timeLoop().timeStepSize()*fluidMatrixInteraction.dsw_dpc(curVolVars.capillaryPressure());
184 }
185
198 template<class PartialDerivativeMatrix>
199 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
200 const Problem& problem,
201 const Element& element,
202 const FVElementGeometry& fvGeometry,
203 const VolumeVariables& curVolVars,
204 const SubControlVolume& scv) const
205 { /* TODO maybe forward to problem for the user to implement the source derivatives?*/ }
206
221 template<class PartialDerivativeMatrices, class T = TypeTag>
222 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::cctpfa, void>
223 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
224 const Problem& problem,
225 const Element& element,
226 const FVElementGeometry& fvGeometry,
227 const ElementVolumeVariables& curElemVolVars,
228 const ElementFluxVariablesCache& elemFluxVarsCache,
229 const SubControlVolumeFace& scvf) const
230 {
231 static_assert(!FluidSystem::isCompressible(0),
232 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
233 static_assert(FluidSystem::viscosityIsConstant(0),
234 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
235
236 // get references to the two participating vol vars & parameters
237 const auto insideScvIdx = scvf.insideScvIdx();
238 const auto outsideScvIdx = scvf.outsideScvIdx();
239 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
240 const auto& insideScv = fvGeometry.scv(insideScvIdx);
241 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
242 const auto& insideVolVars = curElemVolVars[insideScvIdx];
243 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
244
245 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
246 static const auto rho = insideVolVars.density(0);
247 static const auto mu = insideVolVars.viscosity(0);
248 static const auto rho_mu = rho/mu;
249
250 // upwind term
251 // evaluate the current wetting phase Darcy flux and resulting upwind weights
253 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
254 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
255 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
256 const auto outsideWeight = 1.0 - insideWeight;
257 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
258
259 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
260 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(outsideElement, outsideScv, InvalidElemSol{});
261
262 // material law derivatives
263 const auto insideSw = insideVolVars.saturation(0);
264 const auto outsideSw = outsideVolVars.saturation(0);
265 const auto insidePc = insideVolVars.capillaryPressure();
266 const auto outsidePc = outsideVolVars.capillaryPressure();
267 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
268 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
269 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
270 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
271
272 // the transmissibility
273 const auto tij = elemFluxVarsCache[scvf].advectionTij();
274
275 // get references to the two participating derivative matrices
276 auto& dI_dI = derivativeMatrices[insideScvIdx];
277 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
278
279 // partial derivative of the wetting phase flux w.r.t. pw
280 dI_dI[conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
281 dI_dJ[conti0EqIdx][0] += -tij*upwindTerm + rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
282 }
283
295 template<class JacobianMatrix, class T = TypeTag>
296 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
297 addFluxDerivatives(JacobianMatrix& A,
298 const Problem& problem,
299 const Element& element,
300 const FVElementGeometry& fvGeometry,
301 const ElementVolumeVariables& curElemVolVars,
302 const ElementFluxVariablesCache& elemFluxVarsCache,
303 const SubControlVolumeFace& scvf) const
304 {
305 static_assert(!FluidSystem::isCompressible(0),
306 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
307 static_assert(FluidSystem::viscosityIsConstant(0),
308 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
309
310 // get references to the two participating vol vars & parameters
311 const auto insideScvIdx = scvf.insideScvIdx();
312 const auto outsideScvIdx = scvf.outsideScvIdx();
313 const auto& insideScv = fvGeometry.scv(insideScvIdx);
314 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
315 const auto& insideVolVars = curElemVolVars[insideScvIdx];
316 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
317
318 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
319 static const auto rho = insideVolVars.density(0);
320 static const auto mu = insideVolVars.viscosity(0);
321 static const auto rho_mu = rho/mu;
322
323 // upwind term
324 // evaluate the current wetting phase Darcy flux and resulting upwind weights
326 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
327 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
328 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
329 const auto outsideWeight = 1.0 - insideWeight;
330 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
331
332 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
333 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, outsideScv, InvalidElemSol{});
334
335 // material law derivatives
336 const auto insideSw = insideVolVars.saturation(0);
337 const auto outsideSw = outsideVolVars.saturation(0);
338 const auto insidePc = insideVolVars.capillaryPressure();
339 const auto outsidePc = outsideVolVars.capillaryPressure();
340 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
341 const auto dkrw_dsw_outside = outsideFluidMatrixInteraction.dkrw_dsw(outsideSw);
342 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
343 const auto dsw_dpw_outside = -outsideFluidMatrixInteraction.dsw_dpc(outsidePc);
344
345 // so far it was the same as for tpfa
346 // the transmissibilities (flux derivatives with respect to all pw-dofs on the element)
347 const auto ti = AdvectionType::calculateTransmissibilities(
348 problem, element, fvGeometry, curElemVolVars, scvf, elemFluxVarsCache[scvf]
349 );
350
351 // get the rows of the jacobian matrix for the inside/outside scv
352 auto& dI_dJ_inside = A[insideScv.dofIndex()];
353 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
354
355 // add the partial derivatives w.r.t all scvs in the element
356 for (const auto& scvJ : scvs(fvGeometry))
357 {
358 // the transmissibily associated with the scvJ
359 const auto& tj = ti[scvJ.indexInElement()];
360 const auto globalJ = scvJ.dofIndex();
361
362 // partial derivative of the wetting phase flux w.r.t. p_w
363 dI_dJ_inside[globalJ][conti0EqIdx][0] += tj*upwindTerm;
364 dI_dJ_outside[globalJ][conti0EqIdx][0] += -tj*upwindTerm;
365 }
366
367 const auto upwindContributionInside = rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
368 const auto upwindContributionOutside = rho_mu*flux*outsideWeight*dkrw_dsw_outside*dsw_dpw_outside;
369
370 // additional contribution of the upwind term only for inside and outside dof
371 A[insideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] += upwindContributionInside;
372 A[insideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] += upwindContributionOutside;
373
374 A[outsideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionInside;
375 A[outsideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionOutside;
376 }
377
392 template<class PartialDerivativeMatrices>
393 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
394 const Problem& problem,
395 const Element& element,
396 const FVElementGeometry& fvGeometry,
397 const ElementVolumeVariables& curElemVolVars,
398 const ElementFluxVariablesCache& elemFluxVarsCache,
399 const SubControlVolumeFace& scvf) const
400 {
401 static_assert(!FluidSystem::isCompressible(0),
402 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
403 static_assert(FluidSystem::viscosityIsConstant(0),
404 "richards/localresidual.hh: Analytic Jacobian only supports fluids with constant viscosity!");
405
406
407 // get references to the two participating vol vars & parameters
408 const auto insideScvIdx = scvf.insideScvIdx();
409 const auto& insideScv = fvGeometry.scv(insideScvIdx);
410 const auto& insideVolVars = curElemVolVars[insideScvIdx];
411 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
412 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
413
414 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
415 static const auto rho = insideVolVars.density(0);
416 static const auto mu = insideVolVars.viscosity(0);
417 static const auto rho_mu = rho/mu;
418
419 // upwind term
420 // evaluate the current wetting phase Darcy flux and resulting upwind weights
422 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
423 const auto flux = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
424 const auto insideWeight = std::signbit(flux) ? (1.0 - upwindWeight) : upwindWeight;
425 const auto outsideWeight = 1.0 - insideWeight;
426 const auto upwindTerm = rho*insideVolVars.mobility(0)*insideWeight + rho*outsideVolVars.mobility(0)*outsideWeight;
427
428 // material law derivatives
429 const auto insideSw = insideVolVars.saturation(0);
430 const auto insidePc = insideVolVars.capillaryPressure();
431 const auto dkrw_dsw_inside = insideFluidMatrixInteraction.dkrw_dsw(insideSw);
432 const auto dsw_dpw_inside = -insideFluidMatrixInteraction.dsw_dpc(insidePc);
433
434 // the transmissibility
435 const auto tij = elemFluxVarsCache[scvf].advectionTij();
436
437 // partial derivative of the wetting phase flux w.r.t. pw
438 derivativeMatrices[insideScvIdx][conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
439 }
440
452 template<class PartialDerivativeMatrices>
453 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
454 const Problem& problem,
455 const Element& element,
456 const FVElementGeometry& fvGeometry,
457 const ElementVolumeVariables& curElemVolVars,
458 const ElementFluxVariablesCache& elemFluxVarsCache,
459 const SubControlVolumeFace& scvf) const
460 {
461 if constexpr(Detail::hasAddRobinFluxDerivatives<Problem,
462 PartialDerivativeMatrices&, Element, FVElementGeometry,
463 ElementVolumeVariables, ElementFluxVariablesCache, SubControlVolumeFace>()
464 )
465 problem.addRobinFluxDerivatives(derivativeMatrices, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
466 }
467
468private:
469 Implementation *asImp_()
470 { return static_cast<Implementation *> (this); }
471
472 const Implementation *asImp_() const
473 { return static_cast<const Implementation *> (this); }
474};
475
476} // end namespace Dumux
477
478#endif
Element-wise calculation of the Jacobian matrix for problems using the Richards fully implicit models...
Definition: porousmediumflow/richards/localresidual.hh:46
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:199
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:223
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:393
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:100
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:167
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:129
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:453
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:297
Defines all properties used in Dumux.
Type traits.
Helper classes to compute the integration elements.
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:34
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Distance implementation details.
Definition: cvfelocalresidual.hh:25
decltype(std::declval< T >().addRobinFluxDerivatives(std::declval< Ts >()...)) RobinDerivDetector
Definition: porousmediumflow/richards/localresidual.hh:29
static constexpr bool hasAddRobinFluxDerivatives()
Definition: porousmediumflow/richards/localresidual.hh:32
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
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 reference frameworks and formulations available for splitting total fluxes into a advective and d...
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24