PLaSK library
Loading...
Searching...
No Matches
therm_python.cpp
Go to the documentation of this file.
1
/*
2
* This file is part of PLaSK (https://plask.app) by Photonics Group at TUL
3
* Copyright (c) 2022 Lodz University of Technology
4
*
5
* This program is free software: you can redistribute it and/or modify
6
* it under the terms of the GNU General Public License as published by
7
* the Free Software Foundation, version 3.
8
*
9
* This program is distributed in the hope that it will be useful,
10
* but WITHOUT ANY WARRANTY; without even the implied warranty of
11
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12
* GNU General Public License for more details.
13
*/
14
#include <cmath>
15
#include <plask/python.hpp>
16
#include <plask/common/fem/python.hpp>
17
18
using namespace
plask
;
19
using namespace
plask::python
;
20
21
#include "../therm2d.hpp"
22
#include "../therm3d.hpp"
23
using namespace
plask::thermal::tstatic
;
24
25
namespace
plask
{
namespace
thermal {
namespace
tstatic {
26
27
template
<
typename
T>
28
struct
Bc
{
29
30
static
const
char
*
NAME
;
31
static
const
char
*
FIRST
;
32
33
inline
static
double
&
first
(T& self);
34
35
static
void
*
convertible
(
PyObject
* obj) {
36
if
(!
PyDict_Check
(obj))
return
nullptr
;
37
return
obj;
38
}
39
static
void
construct
(
PyObject
* obj, boost::python::converter::rvalue_from_python_stage1_data* data) {
40
void
*
storage
= ((boost::python::converter::rvalue_from_python_storage<Tensor2<double>>*)data)->storage.bytes;
41
double
first
= py::extract<double>(
PyDict_GetItemString
(obj,
FIRST
));
42
PyObject
*
pyambient
=
PyDict_GetItemString
(obj,
"ambient"
);
43
double
ambient =
pyambient
? py::extract<double>(
pyambient
) : 300.;
44
new
(
storage
) T(
first
, ambient);
45
data->convertible =
storage
;
46
}
47
48
Bc
(
const
char
*
doc
) {
49
py::converter::registry::push_back(&
convertible
, &
construct
, boost::python::type_id<T>());
50
py::class_<T>(
NAME
,
doc
, py::init<double,double>())
51
.def(
"__repr__"
, &
Bc<T>::__repr__
)
52
.def(
"__str__"
, &
Bc<T>::__str__
)
53
.def(
"__getitem__"
, &
Bc<T>::__getitem__
)
54
.def(
"__setitem__"
, &
Bc<T>::__setitem__
)
55
;
56
}
57
58
static
std::string
__str__
(T& self) {
59
return
str
(
first
(self)) +
" ("
+
str
(self.ambient) +
"K)"
;
60
}
61
62
static
std::string
__repr__
(T& self) {
63
return
"{'"
+ std::string(
FIRST
) +
"': "
+
str
(
first
(self)) +
", 'ambient': "
+
str
(self.ambient) +
"}"
;
64
}
65
66
static
double
__getitem__
(T& self,
const
std::string& key) {
67
if
(key ==
FIRST
)
return
first
(self);
68
else
if
(key ==
"ambient"
)
return
self.ambient;
69
else
throw
KeyError
(key);
70
}
71
72
static
void
__setitem__
(T& self,
const
std::string& key,
double
value) {
73
if
(key ==
FIRST
)
first
(self) = value;
74
else
if
(key ==
"ambient"
) self.ambient = value;
75
else
throw
KeyError
(key);
76
}
77
78
};
79
80
template
<>
const
char
*
Bc<Convection>::NAME
=
"Convection"
;
81
template
<>
const
char
*
Bc<Convection>::FIRST
=
"coeff"
;
82
template
<>
double
&
Bc<Convection>::first
(
Convection
& self) {
return
self.
coeff
; }
83
84
template
<>
const
char
*
Bc<Radiation>::NAME
=
"Radiation"
;
85
template
<>
const
char
*
Bc<Radiation>::FIRST
=
"emissivity"
;
86
template
<>
double
&
Bc<Radiation>::first
(
Radiation
& self) {
return
self.
emissivity
; }
87
88
}}}
89
90
97
BOOST_PYTHON_MODULE
(
static
)
98
{
99
Bc<Convection>
(
u8
"Convective boundary condition value."
);
100
Bc<Radiation>
(
u8
"Radiative boundary condition value."
);
101
102
{
CLASS
(
ThermalFem2DSolver<Geometry2DCartesian>
,
"Static2D"
,
103
u8
"Finite element thermal solver for 2D Cartesian Geometry."
)
104
METHOD
(compute, compute,
u8
"Run thermal calculations"
, py::arg(
"loops"
)=0);
105
RO_PROPERTY
(err, getErr,
u8
"Maximum estimated error"
);
106
RECEIVER
(inHeat,
""
);
107
PROVIDER
(outTemperature,
""
);
108
PROVIDER
(outHeatFlux,
""
);
109
PROVIDER
(outThermalConductivity,
""
);
110
BOUNDARY_CONDITIONS
(temperature_boundary,
u8
"Boundary conditions for the constant temperature"
);
111
BOUNDARY_CONDITIONS
(heatflux_boundary,
u8
"Boundary conditions for the constant heat flux"
);
112
BOUNDARY_CONDITIONS
(convection_boundary,
u8
"Convective boundary conditions"
);
113
BOUNDARY_CONDITIONS
(radiation_boundary,
u8
"Radiative boundary conditions"
);
114
RW_FIELD
(inittemp,
u8
"Initial temperature"
);
115
RW_FIELD
(maxerr,
u8
"Limit for the temperature updates"
);
116
registerFemSolverWithMaskedMesh
(solver);
117
}
118
119
{
CLASS
(
ThermalFem2DSolver<Geometry2DCylindrical>
,
"StaticCyl"
,
120
u8
"Finite element thermal solver for 2D cylindrical Geometry."
)
121
METHOD
(compute, compute,
u8
"Run thermal calculations"
, py::arg(
"loops"
)=0);
122
RO_PROPERTY
(err, getErr,
u8
"Maximum estimated error"
);
123
RECEIVER
(inHeat,
""
);
124
PROVIDER
(outTemperature,
""
);
125
PROVIDER
(outHeatFlux,
""
);
126
PROVIDER
(outThermalConductivity,
""
);
127
BOUNDARY_CONDITIONS
(temperature_boundary,
u8
"Boundary conditions for the constant temperature"
);
128
BOUNDARY_CONDITIONS
(heatflux_boundary,
u8
"Boundary conditions for the constant heat flux"
);
129
BOUNDARY_CONDITIONS
(convection_boundary,
u8
"Convective boundary conditions"
);
130
BOUNDARY_CONDITIONS
(radiation_boundary,
u8
"Radiative boundary conditions"
);
131
RW_FIELD
(inittemp,
u8
"Initial temperature"
);
132
RW_FIELD
(maxerr,
u8
"Limit for the temperature updates"
);
133
registerFemSolverWithMaskedMesh
(solver);
134
}
135
136
{
CLASS
(
ThermalFem3DSolver
,
"Static3D"
,
u8
"Finite element thermal solver for 3D Geometry."
)
137
METHOD
(compute, compute,
u8
"Run thermal calculations"
, py::arg(
"loops"
)=0);
138
RO_PROPERTY
(err, getErr,
u8
"Maximum estimated error"
);
139
RECEIVER
(inHeat,
""
);
140
solver.setattr(
"inHeatDensity"
, solver.attr(
"inHeat"
));
141
PROVIDER
(outTemperature,
""
);
142
PROVIDER
(outHeatFlux,
""
);
143
PROVIDER
(outThermalConductivity,
""
);
144
BOUNDARY_CONDITIONS
(temperature_boundary,
u8
"Boundary conditions for the constant temperature"
);
145
BOUNDARY_CONDITIONS
(heatflux_boundary,
u8
"Boundary conditions for the constant heat flux"
);
146
BOUNDARY_CONDITIONS
(convection_boundary,
u8
"Convective boundary conditions"
);
147
BOUNDARY_CONDITIONS
(radiation_boundary,
u8
"Radiative boundary conditions"
);
148
RW_FIELD
(inittemp,
u8
"Initial temperature"
);
149
RW_FIELD
(maxerr,
u8
"Limit for the temperature updates"
);
150
registerFemSolverWithMaskedMesh
(solver);
151
}
152
}
solvers
thermal
static
python
therm_python.cpp
Generated by
1.9.8