version 3.8
python/common/volumevariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
13#define DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
14
15#include <dune/common/std/type_traits.hh>
16#include <dune/python/pybind11/pybind11.h>
17#include <dune/python/pybind11/stl.h>
18
19#include <dune/python/common/typeregistry.hh>
20#include <dune/common/classname.hh>
21
22namespace Dumux::Python::Impl {
23
24// helper structs and functions for detecting member functions in volVars
25template <class VolumeVariables>
26using PhaseTemperatureDetector = decltype(std::declval<VolumeVariables>().temperature(0));
27
28template<class VolumeVariables>
29static constexpr bool hasPhaseTemperature()
30{ return Dune::Std::is_detected<PhaseTemperatureDetector, VolumeVariables>::value; }
31
32template <class VolumeVariables>
33using MoleFractionDetector = decltype(std::declval<VolumeVariables>().moleFraction(0, 0));
34
35template<class VolumeVariables>
36static constexpr bool hasMoleFraction()
37{ return Dune::Std::is_detected<MoleFractionDetector, VolumeVariables>::value; }
38
39template <class VolumeVariables>
40using MassFractionDetector = decltype(std::declval<VolumeVariables>().massFraction(0, 0));
41
42template<class VolumeVariables>
43static constexpr bool hasMassFraction()
44{ return Dune::Std::is_detected<MassFractionDetector, VolumeVariables>::value; }
45
46template <class VolumeVariables>
47using SaturationDetector = decltype(std::declval<VolumeVariables>().saturation(0));
48
49template<class VolumeVariables>
50static constexpr bool hasSaturation()
51{ return Dune::Std::is_detected<SaturationDetector, VolumeVariables>::value; }
52
53template <class VolumeVariables>
54using PermeabilityDetector = decltype(std::declval<VolumeVariables>().permeability());
55
56template<class VolumeVariables>
57static constexpr bool hasPermeability()
58{ return Dune::Std::is_detected<PermeabilityDetector, VolumeVariables>::value; }
59
60template<class VolumeVariables>
61void registerVolumeVariables(pybind11::handle scope)
62{
63 using namespace Dune::Python;
64
65 auto [cls, addedToRegistry] = insertClass<VolumeVariables>(
66 scope, "VolumeVariables",
67 GenerateTypeName(Dune::className<VolumeVariables>()),
68 IncludeFiles{}
69 );
70
71 if (addedToRegistry)
72 {
73 using pybind11::operator""_a;
74
75 cls.def("pressure", &VolumeVariables::pressure, "phaseIdx"_a=0);
76 cls.def("density", &VolumeVariables::density, "phaseIdx"_a=0);
77 cls.def("temperature", &VolumeVariables::temperature);
78
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);
89 }
90}
91
92} // namespace Dumux::Python
93
94#endif
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