3.6-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// -*- 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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
25#define DUMUX_PYTHON_COMMON_VOLUMEVARIABLES_HH
26
27#include <dune/common/std/type_traits.hh>
28#include <dune/python/pybind11/pybind11.h>
29#include <dune/python/pybind11/stl.h>
30
31#include <dune/python/common/typeregistry.hh>
32#include <dune/common/classname.hh>
33
34namespace Dumux::Python::Impl {
35
36// helper structs and functions for detecting member functions in volVars
37template <class VolumeVariables>
38using PhaseTemperatureDetector = decltype(std::declval<VolumeVariables>().temperature(0));
39
40template<class VolumeVariables>
41static constexpr bool hasPhaseTemperature()
42{ return Dune::Std::is_detected<PhaseTemperatureDetector, VolumeVariables>::value; }
43
44template <class VolumeVariables>
45using MoleFractionDetector = decltype(std::declval<VolumeVariables>().moleFraction(0, 0));
46
47template<class VolumeVariables>
48static constexpr bool hasMoleFraction()
49{ return Dune::Std::is_detected<MoleFractionDetector, VolumeVariables>::value; }
50
51template <class VolumeVariables>
52using MassFractionDetector = decltype(std::declval<VolumeVariables>().massFraction(0, 0));
53
54template<class VolumeVariables>
55static constexpr bool hasMassFraction()
56{ return Dune::Std::is_detected<MassFractionDetector, VolumeVariables>::value; }
57
58template <class VolumeVariables>
59using SaturationDetector = decltype(std::declval<VolumeVariables>().saturation(0));
60
61template<class VolumeVariables>
62static constexpr bool hasSaturation()
63{ return Dune::Std::is_detected<SaturationDetector, VolumeVariables>::value; }
64
65template <class VolumeVariables>
66using PermeabilityDetector = decltype(std::declval<VolumeVariables>().permeability());
67
68template<class VolumeVariables>
69static constexpr bool hasPermeability()
70{ return Dune::Std::is_detected<PermeabilityDetector, VolumeVariables>::value; }
71
72template<class VolumeVariables>
73void registerVolumeVariables(pybind11::handle scope)
74{
75 using namespace Dune::Python;
76
77 auto [cls, addedToRegistry] = insertClass<VolumeVariables>(
78 scope, "VolumeVariables",
79 GenerateTypeName(Dune::className<VolumeVariables>()),
80 IncludeFiles{}
81 );
82
83 if (addedToRegistry)
84 {
85 using pybind11::operator""_a;
86
87 cls.def("pressure", &VolumeVariables::pressure, "phaseIdx"_a=0);
88 cls.def("density", &VolumeVariables::density, "phaseIdx"_a=0);
89 cls.def("temperature", &VolumeVariables::temperature);
90
91 if constexpr(hasSaturation<VolumeVariables>())
92 cls.def("saturation", &VolumeVariables::saturation, "saturation"_a=0);
93 if constexpr(hasPhaseTemperature<VolumeVariables>())
94 cls.def("temperature", &VolumeVariables::temperature, "phaseIdx"_a=0);
95 if constexpr(hasMoleFraction<VolumeVariables>())
96 cls.def("moleFraction", &VolumeVariables::moleFraction, "phaseIdx"_a, "compIdx"_a);
97 if constexpr(hasMassFraction<VolumeVariables>())
98 cls.def("massFraction", &VolumeVariables::massFraction, "phaseIdx"_a, "compIdx"_a);
99 if constexpr(hasPermeability<VolumeVariables>())
100 cls.def("permeability", &VolumeVariables::permeability);
101 }
102}
103
104} // namespace Dumux::Python
105
106#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