3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
python/common/volumevariables.hh
Go to the documentation of this file.
1#ifndef DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
2#define DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
3
4#include <dune/common/std/type_traits.hh>
5#include <dune/python/pybind11/pybind11.h>
6#include <dune/python/pybind11/stl.h>
7
8#include <dune/python/common/typeregistry.hh>
9#include <dune/common/classname.hh>
10
11namespace Dumux::Python::Impl {
12
13// helper structs and functions for detecting member functions in volVars
14template <class VolumeVariables>
15using PhaseTemperatureDetector = decltype(std::declval<VolumeVariables>().temperature(0));
16
17template<class VolumeVariables>
18static constexpr bool hasPhaseTemperature()
19{ return Dune::Std::is_detected<PhaseTemperatureDetector, VolumeVariables>::value; }
20
21template <class VolumeVariables>
22using MoleFractionDetector = decltype(std::declval<VolumeVariables>().moleFraction(0, 0));
23
24template<class VolumeVariables>
25static constexpr bool hasMoleFraction()
26{ return Dune::Std::is_detected<MoleFractionDetector, VolumeVariables>::value; }
27
28template <class VolumeVariables>
29using MassFractionDetector = decltype(std::declval<VolumeVariables>().massFraction(0, 0));
30
31template<class VolumeVariables>
32static constexpr bool hasMassFraction()
33{ return Dune::Std::is_detected<MassFractionDetector, VolumeVariables>::value; }
34
35template <class VolumeVariables>
36using SaturationDetector = decltype(std::declval<VolumeVariables>().saturation(0));
37
38template<class VolumeVariables>
39static constexpr bool hasSaturation()
40{ return Dune::Std::is_detected<SaturationDetector, VolumeVariables>::value; }
41
42template <class VolumeVariables>
43using PermeabilityDetector = decltype(std::declval<VolumeVariables>().permeability());
44
45template<class VolumeVariables>
46static constexpr bool hasPermeability()
47{ return Dune::Std::is_detected<PermeabilityDetector, VolumeVariables>::value; }
48
49template<class VolumeVariables>
50void registerVolumeVariables(pybind11::handle scope)
51{
52 using namespace Dune::Python;
53
54 auto [cls, addedToRegistry] = insertClass<VolumeVariables>(
55 scope, "VolumeVariables",
56 GenerateTypeName(Dune::className<VolumeVariables>()),
57 IncludeFiles{}
58 );
59
60 if (addedToRegistry)
61 {
62 using pybind11::operator""_a;
63
64 cls.def("pressure", &VolumeVariables::pressure, "phaseIdx"_a=0);
65 cls.def("density", &VolumeVariables::density, "phaseIdx"_a=0);
66 cls.def("temperature", &VolumeVariables::temperature);
67
68 if constexpr(hasSaturation<VolumeVariables>())
69 cls.def("saturation", &VolumeVariables::saturation, "saturation"_a=0);
70 if constexpr(hasPhaseTemperature<VolumeVariables>())
71 cls.def("temperature", &VolumeVariables::temperature, "phaseIdx"_a=0);
72 if constexpr(hasMoleFraction<VolumeVariables>())
73 cls.def("moleFraction", &VolumeVariables::moleFraction, "phaseIdx"_a, "compIdx"_a);
74 if constexpr(hasMassFraction<VolumeVariables>())
75 cls.def("massFraction", &VolumeVariables::massFraction, "phaseIdx"_a, "compIdx"_a);
76 if constexpr(hasPermeability<VolumeVariables>())
77 cls.def("permeability", &VolumeVariables::permeability);
78 }
79}
80
81} // namespace Dumux::Python
82
83#endif
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