14#ifndef DUMUX_RICHARDS_LOCAL_RESIDUAL_HH
15#define DUMUX_RICHARDS_LOCAL_RESIDUAL_HH
28template <
typename T,
typename ...Ts>
29using RobinDerivDetector =
decltype(std::declval<T>().addRobinFluxDerivatives(std::declval<Ts>()...));
31template<
class T,
typename ...Args>
44template<
class TypeTag>
58 using FVElementGeometry =
typename GridGeometry::LocalView;
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;
68 enum { conti0EqIdx = Indices::conti0EqIdx };
74 static constexpr auto liquidPhaseIdx = FluidSystem::phase0Idx;
75 static constexpr auto gasPhaseIdx = FluidSystem::phase1Idx;
76 static constexpr auto liquidCompIdx = FluidSystem::comp0Idx;
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; }
87 using ParentType::ParentType;
101 const SubControlVolume& scv,
102 const VolumeVariables& volVars)
const
105 NumEqVector storage(0.0);
106 storage[conti0EqIdx] = volVars.porosity()
107 * volVars.density(liquidPhaseIdx)
108 * volVars.saturation(liquidPhaseIdx);
111 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, liquidPhaseIdx);
112 EnergyLocalResidual::fluidPhaseStorage(storage, problem, scv, volVars, gasPhaseIdx);
113 EnergyLocalResidual::solidPhaseStorage(storage, scv, volVars);
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const ElementFluxVariablesCache& elemFluxVarsCache)
const
136 FluxVariables fluxVars;
137 fluxVars.init(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
139 NumEqVector flux(0.0);
141 auto upwindTerm = [](
const auto& volVars)
142 {
return volVars.density(liquidPhaseIdx)*volVars.mobility(liquidPhaseIdx); };
144 flux[conti0EqIdx] = fluxVars.advectiveFlux(liquidPhaseIdx, upwindTerm);
147 EnergyLocalResidual::heatConvectionFlux(flux, fluxVars, liquidPhaseIdx);
151 EnergyLocalResidual::heatConductionFlux(flux, fluxVars);
166 template<
class PartialDerivativeMatrix>
168 const Problem& problem,
169 const Element& element,
170 const FVElementGeometry& fvGeometry,
171 const VolumeVariables& curVolVars,
172 const SubControlVolume& scv)
const
174 static_assert(!FluidSystem::isCompressible(0),
175 "richards/localresidual.hh: Analytic Jacobian only supports incompressible fluids!");
177 const auto poreVolume =
Extrusion::volume(fvGeometry, scv)*curVolVars.porosity()*curVolVars.extrusionFactor();
178 static const auto rho = curVolVars.density(0);
182 const auto fluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, scv, InvalidElemSol{});
183 partialDerivatives[conti0EqIdx][0] += -rho*poreVolume/this->timeLoop().timeStepSize()*fluidMatrixInteraction.dsw_dpc(curVolVars.capillaryPressure());
198 template<
class PartialDerivativeMatrix>
200 const Problem& problem,
201 const Element& element,
202 const FVElementGeometry& fvGeometry,
203 const VolumeVariables& curVolVars,
204 const SubControlVolume& scv)
const
221 template<
class PartialDerivativeMatrices,
class T = TypeTag>
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
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!");
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];
246 static const auto rho = insideVolVars.density(0);
247 static const auto mu = insideVolVars.viscosity(0);
248 static const auto rho_mu = rho/mu;
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;
259 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
260 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(outsideElement, outsideScv, InvalidElemSol{});
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);
273 const auto tij = elemFluxVarsCache[scvf].advectionTij();
276 auto& dI_dI = derivativeMatrices[insideScvIdx];
277 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
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;
295 template<
class JacobianMatrix,
class T = TypeTag>
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
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!");
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];
319 static const auto rho = insideVolVars.density(0);
320 static const auto mu = insideVolVars.viscosity(0);
321 static const auto rho_mu = rho/mu;
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;
332 const auto insideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, insideScv, InvalidElemSol{});
333 const auto outsideFluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(element, outsideScv, InvalidElemSol{});
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);
347 const auto ti = AdvectionType::calculateTransmissibilities(
348 problem, element, fvGeometry, curElemVolVars, scvf, elemFluxVarsCache[scvf]
352 auto& dI_dJ_inside = A[insideScv.dofIndex()];
353 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
356 for (
const auto& scvJ : scvs(fvGeometry))
359 const auto& tj = ti[scvJ.indexInElement()];
360 const auto globalJ = scvJ.dofIndex();
363 dI_dJ_inside[globalJ][conti0EqIdx][0] += tj*upwindTerm;
364 dI_dJ_outside[globalJ][conti0EqIdx][0] += -tj*upwindTerm;
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;
371 A[insideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] += upwindContributionInside;
372 A[insideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] += upwindContributionOutside;
374 A[outsideScv.dofIndex()][insideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionInside;
375 A[outsideScv.dofIndex()][outsideScv.dofIndex()][conti0EqIdx][0] -= upwindContributionOutside;
392 template<
class PartialDerivativeMatrices>
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
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!");
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{});
415 static const auto rho = insideVolVars.density(0);
416 static const auto mu = insideVolVars.viscosity(0);
417 static const auto rho_mu = rho/mu;
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;
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);
435 const auto tij = elemFluxVarsCache[scvf].advectionTij();
438 derivativeMatrices[insideScvIdx][conti0EqIdx][0] += tij*upwindTerm + rho_mu*flux*insideWeight*dkrw_dsw_inside*dsw_dpw_inside;
452 template<
class PartialDerivativeMatrices>
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
462 PartialDerivativeMatrices&, Element, FVElementGeometry,
463 ElementVolumeVariables, ElementFluxVariablesCache, SubControlVolumeFace>()
465 problem.addRobinFluxDerivatives(derivativeMatrices, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
469 Implementation *asImp_()
470 {
return static_cast<Implementation *
> (
this); }
472 const Implementation *asImp_()
const
473 {
return static_cast<const Implementation *
> (
this); }
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.
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
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.
Template which always yields a false value.
Definition: common/typetraits/typetraits.hh:24