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;
79 if constexpr(hasSaturation<VolumeVariables>())
81 if constexpr(hasPhaseTemperature<VolumeVariables>())
83 if constexpr(hasMoleFraction<VolumeVariables>())
85 if constexpr(hasMassFraction<VolumeVariables>())
87 if constexpr(hasPermeability<VolumeVariables>())
decltype(std::declval< T >().spatialParams().saturation(std::declval< Ts >()...)) SaturationDetector
Definition: porousmediumflow/tracer/volumevariables.hh:30
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
std::string saturation(int phaseIdx) noexcept
I/O name of saturation for multiphase systems.
Definition: name.hh:31
std::string moleFraction(int phaseIdx, int compIdx) noexcept
I/O name of mole fraction.
Definition: name.hh:98
std::string permeability() noexcept
I/O name of permeability.
Definition: name.hh:131
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:22
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
std::string massFraction(int phaseIdx, int compIdx) noexcept
I/O name of mass fraction.
Definition: name.hh:103