12#ifndef DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
13#define DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
15#include <dune/common/std/type_traits.hh>
16#include <dune/python/pybind11/pybind11.h>
17#include <dune/python/pybind11/stl.h>
19#include <dune/python/common/typeregistry.hh>
20#include <dune/common/classname.hh>
22namespace Dumux::Python::Impl {
25template <
class VolumeVariables>
26using PhaseTemperatureDetector =
decltype(std::declval<VolumeVariables>().temperature(0));
28template<
class VolumeVariables>
29static constexpr bool hasPhaseTemperature()
30{
return Dune::Std::is_detected<PhaseTemperatureDetector, VolumeVariables>::value; }
32template <
class VolumeVariables>
33using MoleFractionDetector =
decltype(std::declval<VolumeVariables>().moleFraction(0, 0));
35template<
class VolumeVariables>
36static constexpr bool hasMoleFraction()
37{
return Dune::Std::is_detected<MoleFractionDetector, VolumeVariables>::value; }
39template <
class VolumeVariables>
40using MassFractionDetector =
decltype(std::declval<VolumeVariables>().massFraction(0, 0));
42template<
class VolumeVariables>
43static constexpr bool hasMassFraction()
44{
return Dune::Std::is_detected<MassFractionDetector, VolumeVariables>::value; }
46template <
class VolumeVariables>
49template<
class VolumeVariables>
50static constexpr bool hasSaturation()
51{
return Dune::Std::is_detected<SaturationDetector, VolumeVariables>::value; }
53template <
class VolumeVariables>
54using PermeabilityDetector =
decltype(std::declval<VolumeVariables>().permeability());
56template<
class VolumeVariables>
57static constexpr bool hasPermeability()
58{
return Dune::Std::is_detected<PermeabilityDetector, VolumeVariables>::value; }
60template<
class VolumeVariables>
61void registerVolumeVariables(pybind11::handle scope)
63 using namespace Dune::Python;
65 auto [cls, addedToRegistry] = insertClass<VolumeVariables>(
66 scope,
"VolumeVariables",
67 GenerateTypeName(Dune::className<VolumeVariables>()),
73 using pybind11::operator
""_a;
75 cls.def(
"pressure", &VolumeVariables::pressure,
"phaseIdx"_a=0);
76 cls.def(
"density", &VolumeVariables::density,
"phaseIdx"_a=0);
77 cls.def(
"temperature", &VolumeVariables::temperature);
79 if constexpr(hasSaturation<VolumeVariables>())
80 cls.def(
"saturation", &VolumeVariables::saturation,
"saturation"_a=0);
81 if constexpr(hasPhaseTemperature<VolumeVariables>())
82 cls.def(
"temperature", &VolumeVariables::temperature,
"phaseIdx"_a=0);
83 if constexpr(hasMoleFraction<VolumeVariables>())
84 cls.def(
"moleFraction", &VolumeVariables::moleFraction,
"phaseIdx"_a,
"compIdx"_a);
85 if constexpr(hasMassFraction<VolumeVariables>())
86 cls.def(
"massFraction", &VolumeVariables::massFraction,
"phaseIdx"_a,
"compIdx"_a);
87 if constexpr(hasPermeability<VolumeVariables>())
88 cls.def(
"permeability", &VolumeVariables::permeability);
decltype(std::declval< T >().spatialParams().saturation(std::declval< Ts >()...)) SaturationDetector
Definition porousmediumflow/tracer/volumevariables.hh:30