3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
deprecated.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 *****************************************************************************/
25#ifndef DUMUX_COMMON_DEPRECATED_HH
26#define DUMUX_COMMON_DEPRECATED_HH
27
28#include <dune/common/deprecated.hh>
32
33namespace Dumux {
34
35#ifndef DOXYGEN // hide from doxygen
36// Helper classes/functions for deprecation
37// Each implementation has to state after which release
38// it will be removed. Implementations in the Deprecated
39// namespace will be removed without
40// deprecation after their usage in the code exprired,
41// so most likely you don't want to use this in your code
42namespace Deprecated {
43
45// Remove the following after Release 3.2 //
47
49// Deprecation warnings for the new material law //
51
52// support old interface of the effective thermal conductivity laws
53template<class E, class SCV, class Sol>
54struct HasNewFIAIF
55{
56 template<class S>
57 auto operator()(S&& sp)
58 -> decltype(sp.fluidMatrixInteraction(std::declval<const E&>(),
59 std::declval<const SCV&>(),
60 std::declval<const Sol&>())) {}
61};
62
63template<class Pos>
64struct HasNewFIAIFAtPos
65{
66 template<class S>
67 auto operator()(S&& sp)
68 -> decltype(sp.fluidMatrixInteractionAtPos(std::declval<const Pos&>())) {}
69};
70
71
72template<class ScalarT, class SpatialParams, class Element, class Scv, class ElemSol>
73class PcKrSwHelper : public FluidMatrix::Adapter<PcKrSwHelper<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::PcKrSw>
74{
75public:
76 using Scalar = ScalarT;
77
78 // pass scalar so template arguments can all be deduced
79 PcKrSwHelper(const Scalar& scalar,
80 const SpatialParams& sp,
81 const Element& element,
82 const Scv& scv,
83 const ElemSol& elemSol)
84 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
85 {}
86
87 Scalar krw(const Scalar sw) const
88 {
89 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
90 return SpatialParams::MaterialLaw::krw(params, sw);
91 }
92
93 Scalar krn(const Scalar sw) const
94 {
95 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
96 return SpatialParams::MaterialLaw::krn(params, sw);
97 }
98
99 Scalar pc(const Scalar sw) const
100 {
101 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
102 return SpatialParams::MaterialLaw::pc(params, sw);
103 }
104
105 Scalar dpc_dsw(const Scalar sw) const
106 {
107 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
108 return SpatialParams::MaterialLaw::dpc_dsw(params, sw);
109 }
110
111 Scalar endPointPc() const
112 {
113 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
114 return SpatialParams::MaterialLaw::endPointPc(params);
115 }
116
117 Scalar sw(const Scalar pc) const
118 {
119 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
120 return SpatialParams::MaterialLaw::sw(params, pc);
121 }
122
123 Scalar dsw_dpc(const Scalar pc) const
124 {
125 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
126 return SpatialParams::MaterialLaw::dsw_dpc(params, pc);
127 }
128
129 Scalar dkrw_dsw(const Scalar sw) const
130 {
131 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
132 return SpatialParams::MaterialLaw::dkrw_dsw(params, sw);
133 }
134
135 Scalar dkrn_dsw(const Scalar sw) const
136 {
137 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
138 return SpatialParams::MaterialLaw::dkrn_dsw(params, sw);
139 }
140
141 const auto& basicParams() const
142 { return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
143
144 const auto& effToAbsParams() const
145 { return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
146
147private:
148 const SpatialParams& spatialParams_;
149 const Element& element_;
150 const Scv& scv_;
151 const ElemSol& elemSol_;
152};
153
154template<class ScalarT, class SpatialParams, class Element, class Scv, class ElemSol>
155class PcKrSwThreePHelper : public FluidMatrix::Adapter<PcKrSwThreePHelper<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::ThreePhasePcKrSw>
156{
157public:
158 using Scalar = ScalarT;
159
160 // pass scalar so template arguments can all be deduced
161 PcKrSwThreePHelper(const Scalar& scalar,
162 const SpatialParams& sp,
163 const Element& element,
164 const Scv& scv,
165 const ElemSol& elemSol)
166 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
167 {}
168
169 Scalar pcgw(const Scalar sw, const Scalar /*dummySn*/) const
170 {
171 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
172 return SpatialParams::MaterialLaw::pcgw(params, sw);
173 }
174
175 Scalar pcnw(const Scalar sw, const Scalar /*dummySn*/) const
176 {
177 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
178 return SpatialParams::MaterialLaw::pcnw(params, sw);
179 }
180
181 Scalar pcgn(const Scalar sw, const Scalar sn) const
182 {
183 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
184 return SpatialParams::MaterialLaw::pcgn(params, sw + sn);
185 }
186
187 Scalar pcAlpha(const Scalar /*dummySw*/, const Scalar sn) const
188 {
189 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
190 return SpatialParams::MaterialLaw::pcAlpha(params, sn);
191 }
192
193 Scalar krw(const Scalar sw, const Scalar sn) const
194 {
195 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
196 return SpatialParams::MaterialLaw::krw(params, sw, sn);
197 }
198
199 Scalar krn(const Scalar sw, const Scalar sn) const
200 {
201 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
202 return SpatialParams::MaterialLaw::krw(params, sw, sn);
203 }
204
205 Scalar krg(const Scalar sw, const Scalar sn) const
206 {
207 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
208 return SpatialParams::MaterialLaw::krg(params, sw, sn);
209 }
210
211 Scalar kr(const int phaseIdx, const Scalar sw, const Scalar sn) const
212 {
213 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
214 return SpatialParams::MaterialLaw::kr(params, phaseIdx, sw, sn, 1 - sw - sn);
215 }
216
217 const auto& basicParams() const
218 { return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
219
220 const auto& effToAbsParams() const
221 { return spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_); }
222
223private:
224 const SpatialParams& spatialParams_;
225 const Element& element_;
226 const Scv& scv_;
227 const ElemSol& elemSol_;
228};
229
230// for implicit models
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,
235 const Scv& scv,
236 const ElemSol& elemSol)
237{
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());
245 else
246 {
247 if constexpr (numPhases == 2)
248 return makeFluidMatrixInteraction(PcKrSwHelper(scalar, sp, element, scv, elemSol));
249 else
250 return makeFluidMatrixInteraction(PcKrSwThreePHelper(scalar, sp, element, scv, elemSol));
251 }
252}
253
254// for sequential models
255template<class Scalar, class SpatialParams, class Element>
256auto makePcKrSw(const Scalar& scalar,
257 const SpatialParams& sp,
258 const Element& element)
259{
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());
264 else
265 {
266 using DummyScv = int;
267 using DummyElemSol = int;
268 return makeFluidMatrixInteraction(PcKrSwHelper(scalar, sp, element, DummyScv(), DummyElemSol()));
269 }
270}
271
272
274// Deprecation warnings for the mp material law stuff //
276
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>
279{
280 using MaterialLaw = typename SpatialParams::MaterialLaw;
281 using MPAdapter = Dumux::MPAdapter<MaterialLaw, NumPhases{}()>;
282public:
283 using Scalar = ScalarT;
284
285 // pass scalar so template arguments can all be deduced
286 PcKrSwMPHelper(const Scalar& scalar,
287 const SpatialParams& sp,
288 const Element& element,
289 const Scv& scv,
290 const ElemSol& elemSol,
291 const NumPhases& np)
292 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
293 {}
294
295 template<class FluidState>
296 auto capillaryPressures(const FluidState& fs, int wPhaseIdx) const
297 {
298 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
299 Dune::FieldVector<Scalar, NumPhases{}()> pc;
300 MPAdapter::capillaryPressures(pc, params, fs, wPhaseIdx);
301 return pc;
302 }
303
304 template<class FluidState>
305 auto relativePermeabilities(const FluidState& fs, int wPhaseIdx) const
306 {
307 const auto& params = spatialParams_.materialLawParamsDeprecated(element_, scv_, elemSol_);
308 Dune::FieldVector<Scalar, NumPhases{}()> kr;
309 MPAdapter::capillaryPressures(kr, params, fs, wPhaseIdx);
310 return kr;
311 }
312
313private:
314 const SpatialParams& spatialParams_;
315 const Element& element_;
316 const Scv& scv_;
317 const ElemSol& elemSol_;
318};
319
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,
324 const Scv& scv,
325 const ElemSol& elemSol,
326 const NumPhases& np)
327{
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());
335 else
336 return makeFluidMatrixInteraction(PcKrSwMPHelper(scalar, sp, element, scv, elemSol, np));
337}
338
340// Deprecation warnings for the kinetic surface areas //
342
343// support old interface of surface area
344template<class E, class SCV, class Sol>
345struct HasNewAns
346{
347 template<class S>
348 auto operator()(S&& sp)
349 -> decltype(sp.nonwettingSolidInterfacialArea(std::declval<const E&>(),
350 std::declval<const SCV&>(),
351 std::declval<const Sol&>())) {}
352};
353
354template<class Pos>
355struct HasNewAnsAtPos
356{
357 template<class S>
358 auto operator()(S&& sp)
359 -> decltype(sp.nonwettingSolidInterfacialAreaAtPos(std::declval<const Pos&>())) {}
360};
361
362// support old interface of surface area
363template<class E, class SCV, class Sol>
364struct HasNewAnw
365{
366 template<class S>
367 auto operator()(S&& sp)
368 -> decltype(sp.wettingNonwettingInterfacialArea(std::declval<const E&>(),
369 std::declval<const SCV&>(),
370 std::declval<const Sol&>())) {}
371};
372
373template<class Pos>
374struct HasNewAnwAtPos
375{
376 template<class S>
377 auto operator()(S&& sp)
378 -> decltype(sp.wettingNonwettingInterfacialAreaAtPos(std::declval<const Pos&>())) {}
379};
380
381// support old interface of surface area
382template<class E, class SCV, class Sol>
383struct HasNewAws
384{
385 template<class S>
386 auto operator()(S&& sp)
387 -> decltype(sp.wettingSolidInterfacialArea(std::declval<const E&>(),
388 std::declval<const SCV&>(),
389 std::declval<const Sol&>())) {}
390};
391
392template<class Pos>
393struct HasNewAwsAtPos
394{
395 template<class S>
396 auto operator()(S&& sp)
397 -> decltype(sp.wettingSolidInterfacialAreaAtPos(std::declval<const Pos&>())) {}
398};
399
400template<class ScalarT, class SpatialParams, class Element, class Scv, class ElemSol>
401class WettingNonwettingInterfacialArea : public FluidMatrix::Adapter<WettingNonwettingInterfacialArea<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::WettingNonwettingInterfacialAreaPcSw>
402{
403public:
404 using Scalar = ScalarT;
405
406 WettingNonwettingInterfacialArea(const Scalar& scalar,
407 const SpatialParams& sp,
408 const Element& element,
409 const Scv& scv,
410 const ElemSol& elemSol)
411 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
412 {}
413
414 const auto& basicParams() const
415 {
416 return spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
417 }
418
419 Scalar area(const Scalar sw, const Scalar pc) const
420 {
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);
425 }
426
427 Scalar darea_dpc(const Scalar sw, const Scalar pc)
428 {
429 const auto& surfaceParams = spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
430 using AwnSurface = typename SpatialParams::AwnSurface;
431 return AwnSurface::dawn_dpc(surfaceParams, sw, pc);
432 }
433
434 Scalar darea_dsw(const Scalar sw, const Scalar pc)
435 {
436 const auto& surfaceParams = spatialParams_.aWettingNonWettingSurfaceParams(element_, scv_, elemSol_);
437 using AwnSurface = typename SpatialParams::AwnSurface;
438 return AwnSurface::dawn_dsw(surfaceParams, sw, pc);
439 }
440
441private:
442 const SpatialParams& spatialParams_;
443 const Element& element_;
444 const Scv& scv_;
445 const ElemSol& elemSol_;
446};
447
448template<class ScalarT, class SpatialParams, class Element, class Scv, class ElemSol>
449class NonwettingSolidInterfacialArea : public FluidMatrix::Adapter<NonwettingSolidInterfacialArea<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::NonwettingSolidInterfacialAreaPcSw>
450{
451public:
452 using Scalar = ScalarT;
453
454 NonwettingSolidInterfacialArea(const Scalar& scalar,
455 const SpatialParams& sp,
456 const Element& element,
457 const Scv& scv,
458 const ElemSol& elemSol)
459 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
460 {}
461
462 const auto& basicParams() const
463 {
464 return spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
465 }
466
467 Scalar area(const Scalar sw, const Scalar pc) const
468 {
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);
473 }
474
475 Scalar darea_dpc(const Scalar sw, const Scalar pc)
476 {
477 const auto& surfaceParams = spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
478 using AnsSurface = typename SpatialParams::AnsSurface;
479 return AnsSurface::dawn_dpc(surfaceParams, sw, pc);
480 }
481
482 Scalar darea_dsw(const Scalar sw, const Scalar pc)
483 {
484 const auto& surfaceParams = spatialParams_.aNonWettingSolidSurfaceParams(element_, scv_, elemSol_);
485 using AnsSurface = typename SpatialParams::AnsSurface;
486 return AnsSurface::dawn_dsw(surfaceParams, sw, pc);
487 }
488
489private:
490 const SpatialParams& spatialParams_;
491 const Element& element_;
492 const Scv& scv_;
493 const ElemSol& elemSol_;
494};
495
496template<class ScalarT, class SpatialParams, class Element, class Scv, class ElemSol>
497class WettingSolidInterfacialArea : public FluidMatrix::Adapter<WettingSolidInterfacialArea<ScalarT, SpatialParams, Element, Scv, ElemSol>, FluidMatrix::WettingSolidInterfacialAreaPcSw>
498{
499public:
500 using Scalar = ScalarT;
501
502 WettingSolidInterfacialArea(const Scalar& scalar,
503 const SpatialParams& sp,
504 const Element& element,
505 const Scv& scv,
506 const ElemSol& elemSol)
507 : spatialParams_(sp), element_(element), scv_(scv), elemSol_(elemSol)
508 {}
509
510 const auto& basicParams() const
511 {
512 return spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
513 }
514
515 Scalar area(const Scalar sw, const Scalar pc) const
516 {
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);
521 }
522
523 Scalar darea_dpc(const Scalar sw, const Scalar pc)
524 {
525 const auto& surfaceParams = spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
526 using AwsSurface = typename SpatialParams::AwsSurface;
527 return AwsSurface::dawn_dpc(surfaceParams, sw, pc);
528 }
529
530 Scalar darea_dsw(const Scalar sw, const Scalar pc)
531 {
532 const auto& surfaceParams = spatialParams_.aWettingSolidSurfaceParams(element_, scv_, elemSol_);
533 using AwsSurface = typename SpatialParams::AwsSurface;
534 return AwsSurface::dawn_dsw(surfaceParams, sw, pc);
535 }
536
537private:
538 const SpatialParams& spatialParams_;
539 const Element& element_;
540 const Scv& scv_;
541 const ElemSol& elemSol_;
542};
543
544template<class Scalar, class SpatialParams, class Element, class Scv, class ElemSol>
545auto makeInterfacialArea(const Scalar& scalar,
546 const SpatialParams& sp,
547 const Element& element,
548 const Scv& scv,
549 const ElemSol& elemSol)
550{
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());
558 else
559 return makeFluidMatrixInteraction(WettingNonwettingInterfacialArea(scalar, sp, element, scv, elemSol),
560 NonwettingSolidInterfacialArea(scalar, sp, element, scv, elemSol),
561 WettingSolidInterfacialArea(scalar, sp, element, scv, elemSol));
562}
563
564} // end namespace Deprecated
565#endif
566
567} // end namespace Dumux
568#endif
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
Definition: adapt.hh:29
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