Skip to content
GitLab
Menu
Projects
Groups
Snippets
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
dumux-repositories
dumux
Commits
e2eb7b40
Commit
e2eb7b40
authored
Oct 29, 2020
by
Timo Koch
Committed by
Kilian Weishaupt
Oct 30, 2020
Browse files
[deprecated] Add deprecation helper for mp material laws
parent
a87c1585
Changes
1
Hide whitespace changes
Inline
Side-by-side
dumux/common/deprecated.hh
View file @
e2eb7b40
...
...
@@ -28,6 +28,7 @@
#include <dune/common/deprecated.hh>
#include <dumux/common/typetraits/isvalid.hh>
#include <dumux/material/fluidmatrixinteractions/fluidmatrixinteraction.hh>
#include <dumux/material/fluidmatrixinteractions/mp/mpadapter.hh>
namespace
Dumux
{
...
...
@@ -188,6 +189,71 @@ auto makePcKrSw(const Scalar& scalar,
}
///////////////////////////////////////////////////////////////
// Deprecation warnings for the mp material law stuff //
///////////////////////////////////////////////////////////////
template
<
class
ScalarT
,
class
SpatialParams
,
class
Element
,
class
Scv
,
class
ElemSol
,
class
NumPhases
>
class
PcKrSwMPHelper
:
public
FluidMatrix
::
Adapter
<
PcKrSwMPHelper
<
ScalarT
,
SpatialParams
,
Element
,
Scv
,
ElemSol
,
NumPhases
>
,
FluidMatrix
::
MultiPhasePcKrSw
>
{
using
MaterialLaw
=
typename
SpatialParams
::
MaterialLaw
;
using
MPAdapter
=
Dumux
::
MPAdapter
<
MaterialLaw
,
NumPhases
{}()
>
;
public:
using
Scalar
=
ScalarT
;
// pass scalar so template arguments can all be deduced
PcKrSwMPHelper
(
const
Scalar
&
scalar
,
const
SpatialParams
&
sp
,
const
Element
&
element
,
const
Scv
&
scv
,
const
ElemSol
&
elemSol
,
const
NumPhases
&
np
)
:
spatialParams_
(
sp
),
element_
(
element
),
scv_
(
scv
),
elemSol_
(
elemSol
)
{}
template
<
class
FluidState
>
auto
capillaryPressures
(
const
FluidState
&
fs
,
int
wPhaseIdx
)
const
{
const
auto
&
params
=
spatialParams_
.
materialLawParamsDeprecated
(
element_
,
scv_
,
elemSol_
);
Dune
::
FieldVector
<
Scalar
,
NumPhases
{}()
>
pc
;
MPAdapter
::
capillaryPressures
(
pc
,
params
,
fs
,
wPhaseIdx
);
return
pc
;
}
template
<
class
FluidState
>
auto
relativePermeabilities
(
const
FluidState
&
fs
,
int
wPhaseIdx
)
const
{
const
auto
&
params
=
spatialParams_
.
materialLawParamsDeprecated
(
element_
,
scv_
,
elemSol_
);
Dune
::
FieldVector
<
Scalar
,
NumPhases
{}()
>
kr
;
MPAdapter
::
capillaryPressures
(
kr
,
params
,
fs
,
wPhaseIdx
);
return
kr
;
}
private:
const
SpatialParams
&
spatialParams_
;
const
Element
&
element_
;
const
Scv
&
scv_
;
const
ElemSol
&
elemSol_
;
};
template
<
class
Scalar
,
class
SpatialParams
,
class
Element
,
class
Scv
,
class
ElemSol
,
class
NumPhases
>
auto
makeMPPcKrSw
(
const
Scalar
&
scalar
,
const
SpatialParams
&
sp
,
const
Element
&
element
,
const
Scv
&
scv
,
const
ElemSol
&
elemSol
,
const
NumPhases
&
np
)
{
using
GlobalPosition
=
typename
Element
::
Geometry
::
GlobalCoordinate
;
constexpr
bool
hasNew
=
decltype
(
isValid
(
HasNewFIAIF
<
Element
,
Scv
,
ElemSol
>
()).
template
check
<
SpatialParams
>())
::
value
;
constexpr
bool
hasNewAtPos
=
decltype
(
isValid
(
HasNewFIAIFAtPos
<
GlobalPosition
>
()).
template
check
<
SpatialParams
>())
::
value
;
if
constexpr
(
hasNew
)
return
sp
.
fluidMatrixInteraction
(
element
,
scv
,
elemSol
);
else
if
constexpr
(
hasNewAtPos
)
return
sp
.
fluidMatrixInteractionAtPos
(
scv
.
center
());
else
return
makeFluidMatrixInteraction
(
PcKrSwMPHelper
(
scalar
,
sp
,
element
,
scv
,
elemSol
,
np
));
}
///////////////////////////////////////////////////////////////
// Deprecation warnings for the kinetic surface areas //
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment