25#ifndef DUMUX_COMMON_DEPRECATED_HH
26#define DUMUX_COMMON_DEPRECATED_HH
28#include <dune/common/deprecated.hh>
53template<
class E,
class SCV,
class Sol>
57 auto operator()(S&& sp)
58 ->
decltype(sp.fluidMatrixInteraction(std::declval<const E&>(),
59 std::declval<const SCV&>(),
60 std::declval<const Sol&>())) {}
64struct HasNewFIAIFAtPos
67 auto operator()(S&& sp)
68 ->
decltype(sp.fluidMatrixInteractionAtPos(std::declval<const Pos&>())) {}
72template<
class ScalarT,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
73class PcKrSwHelper :
public FluidMatrix::Adapter<PcKrSwHelper<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::PcKrSw>
76 using Scalar = ScalarT;
79 PcKrSwHelper(
const Scalar& scalar,
80 const SpatialParams& sp,
81 const Element& element,
83 const ElemSol& elemSol)
84 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
87 Scalar krw(
const Scalar sw)
const
89 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
90 return SpatialParams::MaterialLaw::krw(params, sw);
93 Scalar krn(
const Scalar sw)
const
95 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
96 return SpatialParams::MaterialLaw::krn(params, sw);
99 Scalar pc(
const Scalar sw)
const
101 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
102 return SpatialParams::MaterialLaw::pc(params, sw);
105 Scalar dpc_dsw(
const Scalar sw)
const
107 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
108 return SpatialParams::MaterialLaw::dpc_dsw(params, sw);
111 Scalar endPointPc()
const
113 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
114 return SpatialParams::MaterialLaw::endPointPc(params);
117 Scalar sw(
const Scalar pc)
const
119 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
120 return SpatialParams::MaterialLaw::sw(params, pc);
123 Scalar dsw_dpc(
const Scalar pc)
const
125 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
126 return SpatialParams::MaterialLaw::dsw_dpc(params, pc);
129 Scalar dkrw_dsw(
const Scalar sw)
const
131 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
132 return SpatialParams::MaterialLaw::dkrw_dsw(params, sw);
135 Scalar dkrn_dsw(
const Scalar sw)
const
137 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
138 return SpatialParams::MaterialLaw::dkrn_dsw(params, sw);
141 const auto& basicParams()
const
142 {
return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
144 const auto& effToAbsParams()
const
145 {
return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
148 const SpatialParams& spatialParams_;
149 const Element& element_;
151 const ElemSol& elemSol_;
154template<
class ScalarT,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
155class PcKrSwThreePHelper :
public FluidMatrix::Adapter<PcKrSwThreePHelper<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::ThreePhasePcKrSw>
158 using Scalar = ScalarT;
161 PcKrSwThreePHelper(
const Scalar& scalar,
162 const SpatialParams& sp,
163 const Element& element,
165 const ElemSol& elemSol)
166 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
169 Scalar pcgw(
const Scalar sw,
const Scalar )
const
171 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
172 return SpatialParams::MaterialLaw::pcgw(params, sw);
175 Scalar pcnw(
const Scalar sw,
const Scalar )
const
177 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
178 return SpatialParams::MaterialLaw::pcnw(params, sw);
181 Scalar pcgn(
const Scalar sw,
const Scalar sn)
const
183 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
184 return SpatialParams::MaterialLaw::pcgn(params, sw + sn);
187 Scalar pcAlpha(
const Scalar ,
const Scalar sn)
const
189 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
190 return SpatialParams::MaterialLaw::pcAlpha(params, sn);
193 Scalar krw(
const Scalar sw,
const Scalar sn)
const
195 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
196 return SpatialParams::MaterialLaw::krw(params, sw, sn);
199 Scalar krn(
const Scalar sw,
const Scalar sn)
const
201 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
202 return SpatialParams::MaterialLaw::krw(params, sw, sn);
205 Scalar krg(
const Scalar sw,
const Scalar sn)
const
207 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
208 return SpatialParams::MaterialLaw::krg(params, sw, sn);
211 Scalar kr(
const int phaseIdx,
const Scalar sw,
const Scalar sn)
const
213 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
214 return SpatialParams::MaterialLaw::kr(params, phaseIdx, sw, sn, 1 - sw - sn);
217 const auto& basicParams()
const
218 {
return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
220 const auto& effToAbsParams()
const
221 {
return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
224 const SpatialParams& spatialParams_;
225 const Element& element_;
227 const ElemSol& elemSol_;
231template<
int numPhases = 2,
class Scalar,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
232auto makePcKrSw(
const Scalar& scalar,
233 const SpatialParams& sp,
234 const Element& element,
236 const ElemSol& elemSol)
238 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
239 constexpr bool hasNew =
decltype(
isValid(HasNewFIAIF<Element, Scv, ElemSol>()).template check<SpatialParams>())::value;
240 [[maybe_unused]]
constexpr bool hasNewAtPos =
decltype(
isValid(HasNewFIAIFAtPos<GlobalPosition>()).template check<SpatialParams>())::value;
241 if constexpr (hasNew)
242 return sp.fluidMatrixInteraction(element, scv, elemSol);
243 else if constexpr (hasNewAtPos)
244 return sp.fluidMatrixInteractionAtPos(scv.center());
247 if constexpr (numPhases == 2)
255template<
class Scalar,
class SpatialParams,
class Element>
256auto makePcKrSw(
const Scalar& scalar,
257 const SpatialParams& sp,
258 const Element& element)
260 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
261 constexpr bool hasNewAtPos =
decltype(
isValid(HasNewFIAIFAtPos<GlobalPosition>()).template check<SpatialParams>())::value;
262 if constexpr (hasNewAtPos)
263 return sp.fluidMatrixInteractionAtPos(element.geometry().center());
266 using DummyScv = int;
267 using DummyElemSol = int;
277template<
class ScalarT,
class SpatialParams,
class Element,
class Scv,
class ElemSol,
class NumPhases>
278class PcKrSwMPHelper :
public FluidMatrix::Adapter<PcKrSwMPHelper<ScalarT, SpatialParams, Element, Scv, ElemSol, NumPhases>, FluidMatrix::MultiPhasePcKrSw>
280 using MaterialLaw =
typename SpatialParams::MaterialLaw;
283 using Scalar = ScalarT;
286 PcKrSwMPHelper(
const Scalar& scalar,
287 const SpatialParams& sp,
288 const Element& element,
290 const ElemSol& elemSol,
292 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
295 template<
class Flu
idState>
296 auto capillaryPressures(
const FluidState& fs,
int wPhaseIdx)
const
298 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
299 Dune::FieldVector<Scalar, NumPhases{}()> pc;
300 MPAdapter::capillaryPressures(pc, params, fs, wPhaseIdx);
304 template<
class Flu
idState>
305 auto relativePermeabilities(
const FluidState& fs,
int wPhaseIdx)
const
307 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
308 Dune::FieldVector<Scalar, NumPhases{}()> kr;
309 MPAdapter::capillaryPressures(kr, params, fs, wPhaseIdx);
314 const SpatialParams& spatialParams_;
315 const Element& element_;
317 const ElemSol& elemSol_;
320template<
class Scalar,
class SpatialParams,
class Element,
class Scv,
class ElemSol,
class NumPhases>
321auto makeMPPcKrSw(
const Scalar& scalar,
322 const SpatialParams& sp,
323 const Element& element,
325 const ElemSol& elemSol,
328 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
329 constexpr bool hasNew =
decltype(
isValid(HasNewFIAIF<Element, Scv, ElemSol>()).template check<SpatialParams>())::value;
330 constexpr bool hasNewAtPos =
decltype(
isValid(HasNewFIAIFAtPos<GlobalPosition>()).template check<SpatialParams>())::value;
331 if constexpr (hasNew)
332 return sp.fluidMatrixInteraction(element, scv, elemSol);
333 else if constexpr (hasNewAtPos)
334 return sp.fluidMatrixInteractionAtPos(scv.center());
344template<
class E,
class SCV,
class Sol>
348 auto operator()(S&& sp)
349 ->
decltype(sp.nonwettingSolidInterfacialArea(std::declval<const E&>(),
350 std::declval<const SCV&>(),
351 std::declval<const Sol&>())) {}
358 auto operator()(S&& sp)
359 ->
decltype(sp.nonwettingSolidInterfacialAreaAtPos(std::declval<const Pos&>())) {}
363template<
class E,
class SCV,
class Sol>
367 auto operator()(S&& sp)
368 ->
decltype(sp.wettingNonwettingInterfacialArea(std::declval<const E&>(),
369 std::declval<const SCV&>(),
370 std::declval<const Sol&>())) {}
377 auto operator()(S&& sp)
378 ->
decltype(sp.wettingNonwettingInterfacialAreaAtPos(std::declval<const Pos&>())) {}
382template<
class E,
class SCV,
class Sol>
386 auto operator()(S&& sp)
387 ->
decltype(sp.wettingSolidInterfacialArea(std::declval<const E&>(),
388 std::declval<const SCV&>(),
389 std::declval<const Sol&>())) {}
396 auto operator()(S&& sp)
397 ->
decltype(sp.wettingSolidInterfacialAreaAtPos(std::declval<const Pos&>())) {}
400template<
class ScalarT,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
401class WettingNonwettingInterfacialArea :
public FluidMatrix::Adapter<WettingNonwettingInterfacialArea<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::WettingNonwettingInterfacialAreaPcSw>
404 using Scalar = ScalarT;
406 WettingNonwettingInterfacialArea(
const Scalar& scalar,
407 const SpatialParams& sp,
408 const Element& element,
410 const ElemSol& elemSol)
411 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
414 const auto& basicParams()
const
416 return spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
419 Scalar area(
const Scalar sw,
const Scalar pc)
const
421 const auto& surfaceParams = spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
422 const auto& materialParams = spatialParams_.materialLawParams(element_, scv_, elemSol_);
423 using AwnSurface =
typename SpatialParams::AwnSurface;
424 return AwnSurface::interfacialArea(surfaceParams, materialParams, sw, pc);
427 Scalar darea_dpc(
const Scalar sw,
const Scalar pc)
429 const auto& surfaceParams = spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
430 using AwnSurface =
typename SpatialParams::AwnSurface;
431 return AwnSurface::dawn_dpc(surfaceParams, sw, pc);
434 Scalar darea_dsw(
const Scalar sw,
const Scalar pc)
436 const auto& surfaceParams = spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
437 using AwnSurface =
typename SpatialParams::AwnSurface;
438 return AwnSurface::dawn_dsw(surfaceParams, sw, pc);
442 const SpatialParams& spatialParams_;
443 const Element& element_;
445 const ElemSol& elemSol_;
448template<
class ScalarT,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
449class NonwettingSolidInterfacialArea :
public FluidMatrix::Adapter<NonwettingSolidInterfacialArea<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::NonwettingSolidInterfacialAreaPcSw>
452 using Scalar = ScalarT;
454 NonwettingSolidInterfacialArea(
const Scalar& scalar,
455 const SpatialParams& sp,
456 const Element& element,
458 const ElemSol& elemSol)
459 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
462 const auto& basicParams()
const
464 return spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
467 Scalar area(
const Scalar sw,
const Scalar pc)
const
469 const auto& surfaceParams = spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
470 const auto& materialParams = spatialParams_.materialLawParams(element_, scv_, elemSol_);
471 using AnsSurface =
typename SpatialParams::AnsSurface;
472 return AnsSurface::interfacialArea(surfaceParams, materialParams, sw, pc);
475 Scalar darea_dpc(
const Scalar sw,
const Scalar pc)
477 const auto& surfaceParams = spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
478 using AnsSurface =
typename SpatialParams::AnsSurface;
479 return AnsSurface::dawn_dpc(surfaceParams, sw, pc);
482 Scalar darea_dsw(
const Scalar sw,
const Scalar pc)
484 const auto& surfaceParams = spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
485 using AnsSurface =
typename SpatialParams::AnsSurface;
486 return AnsSurface::dawn_dsw(surfaceParams, sw, pc);
490 const SpatialParams& spatialParams_;
491 const Element& element_;
493 const ElemSol& elemSol_;
496template<
class ScalarT,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
497class WettingSolidInterfacialArea :
public FluidMatrix::Adapter<WettingSolidInterfacialArea<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::WettingSolidInterfacialAreaPcSw>
500 using Scalar = ScalarT;
502 WettingSolidInterfacialArea(
const Scalar& scalar,
503 const SpatialParams& sp,
504 const Element& element,
506 const ElemSol& elemSol)
507 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
510 const auto& basicParams()
const
512 return spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
515 Scalar area(
const Scalar sw,
const Scalar pc)
const
517 const auto& surfaceParams = spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
518 const auto& materialParams = spatialParams_.materialLawParams(element_, scv_, elemSol_);
519 using AwsSurface =
typename SpatialParams::AwsSurface;
520 return AwsSurface::interfacialArea(surfaceParams, materialParams, sw, pc);
523 Scalar darea_dpc(
const Scalar sw,
const Scalar pc)
525 const auto& surfaceParams = spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
526 using AwsSurface =
typename SpatialParams::AwsSurface;
527 return AwsSurface::dawn_dpc(surfaceParams, sw, pc);
530 Scalar darea_dsw(
const Scalar sw,
const Scalar pc)
532 const auto& surfaceParams = spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
533 using AwsSurface =
typename SpatialParams::AwsSurface;
534 return AwsSurface::dawn_dsw(surfaceParams, sw, pc);
538 const SpatialParams& spatialParams_;
539 const Element& element_;
541 const ElemSol& elemSol_;
544template<
class Scalar,
class SpatialParams,
class Element,
class Scv,
class ElemSol>
545auto makeInterfacialArea(
const Scalar& scalar,
546 const SpatialParams& sp,
547 const Element& element,
549 const ElemSol& elemSol)
551 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
552 constexpr bool hasNew =
decltype(
isValid(HasNewFIAIF<Element, Scv, ElemSol>()).template check<SpatialParams>())::value;
553 constexpr bool hasNewAtPos =
decltype(
isValid(HasNewFIAIFAtPos<GlobalPosition>()).template check<SpatialParams>())::value;
554 if constexpr (hasNew)
555 return sp.fluidMatrixInteraction(element, scv, elemSol);
556 else if constexpr (hasNewAtPos)
557 return sp.fluidMatrixInteractionAtPos(scv.center());
560 NonwettingSolidInterfacialArea(scalar, sp, element, scv, elemSol),
561 WettingSolidInterfacialArea(scalar, sp, element, scv, elemSol));
A helper function for class member function introspection.
Makes the capillary pressure-saturation relations available under the M-phase API for material laws.
Wrapper type to combine an arbitrary number of different laws for fluid-matrix interaction (e....
auto makeFluidMatrixInteraction(Laws &&... laws)
Helper function to create an FluidMatrixInteraction object containing an arbitrary number of fluid ma...
Definition: fluidmatrixinteraction.hh:51
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
An adapter for mpnc to use the capillary pressure-saturation relationships.
Definition: mpadapter.hh:45