24#ifndef DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
25#define DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
27#include <dune/common/std/type_traits.hh>
28#include <dune/python/pybind11/pybind11.h>
29#include <dune/python/pybind11/stl.h>
31#include <dune/python/common/typeregistry.hh>
32#include <dune/common/classname.hh>
34namespace Dumux::Python::Impl {
37template <
class VolumeVariables>
38using PhaseTemperatureDetector =
decltype(std::declval<VolumeVariables>().temperature(0));
40template<
class VolumeVariables>
41static constexpr bool hasPhaseTemperature()
42{
return Dune::Std::is_detected<PhaseTemperatureDetector, VolumeVariables>::value; }
44template <
class VolumeVariables>
45using MoleFractionDetector =
decltype(std::declval<VolumeVariables>().moleFraction(0, 0));
47template<
class VolumeVariables>
48static constexpr bool hasMoleFraction()
49{
return Dune::Std::is_detected<MoleFractionDetector, VolumeVariables>::value; }
51template <
class VolumeVariables>
52using MassFractionDetector =
decltype(std::declval<VolumeVariables>().massFraction(0, 0));
54template<
class VolumeVariables>
55static constexpr bool hasMassFraction()
56{
return Dune::Std::is_detected<MassFractionDetector, VolumeVariables>::value; }
58template <
class VolumeVariables>
61template<
class VolumeVariables>
62static constexpr bool hasSaturation()
63{
return Dune::Std::is_detected<SaturationDetector, VolumeVariables>::value; }
65template <
class VolumeVariables>
66using PermeabilityDetector =
decltype(std::declval<VolumeVariables>().permeability());
68template<
class VolumeVariables>
69static constexpr bool hasPermeability()
70{
return Dune::Std::is_detected<PermeabilityDetector, VolumeVariables>::value; }
72template<
class VolumeVariables>
73void registerVolumeVariables(pybind11::handle scope)
75 using namespace Dune::Python;
77 auto [cls, addedToRegistry] = insertClass<VolumeVariables>(
78 scope,
"VolumeVariables",
79 GenerateTypeName(Dune::className<VolumeVariables>()),
85 using pybind11::operator
""_a;
91 if constexpr(hasSaturation<VolumeVariables>())
93 if constexpr(hasPhaseTemperature<VolumeVariables>())
95 if constexpr(hasMoleFraction<VolumeVariables>())
97 if constexpr(hasMassFraction<VolumeVariables>())
99 if constexpr(hasPermeability<VolumeVariables>())
decltype(std::declval< T >().spatialParams().saturation(std::declval< Ts >()...)) SaturationDetector
Definition: porousmediumflow/tracer/volumevariables.hh:42
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:43
std::string moleFraction(int phaseIdx, int compIdx) noexcept
I/O name of mole fraction.
Definition: name.hh:110
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:143
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
std::string massFraction(int phaseIdx, int compIdx) noexcept
I/O name of mass fraction.
Definition: name.hh:115