PLaSK library
Loading...
Searching...
No Matches
electr_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>
17using namespace plask;
18using namespace plask::python;
19
20#include "../beta.hpp"
21#include "../electr2d.hpp"
22#include "../electr3d.hpp"
23using namespace plask::electrical::shockley;
24
25static py::object outPotential(const py::object& self) {
26 throw TypeError(u8"{}: 'outPotential' is reserved for drift-diffusion model; use 'outVoltage' instead",
27 std::string(py::extract<std::string>(self.attr("id"))));
28 return py::object();
29}
30
31template <typename GeometryT> struct Shockley : BetaSolver<GeometryT> {
32 std::vector<py::object> beta_function, js_function;
33
34 Shockley(const std::string& id = "") : BetaSolver<GeometryT>(id) {}
35
36 py::object getBeta0() const { return getBeta(0); }
37 void setBeta0(const py::object& value) { setBeta(0, value); }
38
39 py::object getJs0() const { return getJs(0); }
40 void setJs0(const py::object& value) { setJs(0, value); }
41
42 py::object __getattr__(const std::string& attr) const {
43 try {
44 if (attr.substr(0, 4) == "beta") return py::object(getBeta(boost::lexical_cast<size_t>(attr.substr(4))));
45 if (attr.substr(0, 2) == "js") return py::object(getJs(boost::lexical_cast<size_t>(attr.substr(2))));
46 } catch (boost::bad_lexical_cast&) {
47 }
48 throw AttributeError(u8"'{0}' object has no attribute '{1}'", this->getClassName(), attr);
49 }
50
51 static void __setattr__(const py::object& oself, const std::string& attr, const py::object& value) {
52 Shockley<GeometryT>& self = py::extract<Shockley<GeometryT>&>(oself);
53 try {
54 if (attr.substr(0, 4) == "beta") {
55 self.setBeta(boost::lexical_cast<size_t>(attr.substr(4)), value);
56 return;
57 }
58 if (attr.substr(0, 2) == "js") {
59 self.setJs(boost::lexical_cast<size_t>(attr.substr(2)), value);
60 return;
61 }
62 } catch (boost::bad_lexical_cast&) {
63 }
64 oself.attr("__class__").attr("__base__").attr("__setattr__")(oself, attr, value);
65 }
66
67 py::object getBeta(size_t n) const {
68 if (n < beta_function.size() && !beta_function[n].is_none()) return beta_function[n];
69 return py::object(BetaSolver<GeometryT>::getBeta(n));
70 }
71
72 void setBeta(size_t n, const py::object& value) {
73 py::extract<double> val(value);
74 if (val.check()) {
76 } else if (PyCallable_Check(value.ptr())) {
77 if (this->beta_function.size() <= n) this->beta_function.resize(n + 1);
78 beta_function[n] = value;
79 this->invalidate();
80 } else {
81 throw TypeError("{}: beta{} must be a float or a callable", this->getId(), n);
82 }
83 }
84
85 py::object getJs(size_t n) const {
86 if (n < js_function.size() && !js_function[n].is_none()) return js_function[n];
87 return py::object(BetaSolver<GeometryT>::getJs(n));
88 }
89
90 void setJs(size_t n, const py::object& value) {
91 py::extract<double> val(value);
92 if (val.check()) {
94 } else if (PyCallable_Check(value.ptr())) {
95 if (this->js_function.size() <= n) this->js_function.resize(n + 1);
96 js_function[n] = value;
97 this->invalidate();
98 } else {
99 throw TypeError("{}: js{} must be a float or a callable", this->getId(), n);
100 }
101 }
102
103 Tensor2<double> activeCond(size_t n, double U, double jy, double T) override {
104 double beta = (n < beta_function.size() && !beta_function[n].is_none()) ? py::extract<double>(beta_function[n](T))
106 double js = (n < js_function.size() && !js_function[n].is_none()) ? py::extract<double>(js_function[n](T))
108 jy = abs(jy);
109 return Tensor2<double>(0., 10. * jy * beta * this->active[n].height / log(1e7 * jy / js + 1.));
110 }
111};
112
113template <typename GeometryT>
114struct PythonCondSolver : public std::conditional<std::is_same<GeometryT, Geometry3D>::value,
115 ElectricalFem3DSolver,
116 ElectricalFem2DSolver<GeometryT>>::type {
117 typedef typename std::conditional<std::is_same<GeometryT, Geometry3D>::value,
120
121 std::vector<py::object> cond_function;
122
123 PythonCondSolver(const std::string& id = "") : BaseClass(id) {}
124
125 py::object getCond0() const { return getCond(0); }
126 void setCond0(const py::object& value) { setCond(0, value); }
127
128 py::object getCond(size_t n) const {
129 if (n < cond_function.size()) return cond_function[n];
130 return py::object();
131 }
132
133 void setCond(size_t n, const py::object& value) {
134 if (PyCallable_Check(value.ptr())) {
135 if (this->cond_function.size() <= n) this->cond_function.resize(n + 1);
136 cond_function[n] = value;
137 this->invalidate();
138 } else {
139 throw TypeError("{}: cond{} must be a a callable", this->getId(), n);
140 }
141 }
142
143 py::object __getattr__(const std::string& attr) const {
144 try {
145 if (attr.substr(0, 4) == "cond") return py::object(getCond(boost::lexical_cast<size_t>(attr.substr(4))));
146 } catch (boost::bad_lexical_cast&) {
147 }
148 throw AttributeError(u8"'{0}' object has no attribute '{1}'", this->getClassName(), attr);
149 }
150
151 static void __setattr__(const py::object& oself, const std::string& attr, const py::object& value) {
152 PythonCondSolver<GeometryT>& self = py::extract<PythonCondSolver<GeometryT>&>(oself);
153
154 try {
155 if (attr.substr(0, 4) == "cond") {
156 self.setCond(boost::lexical_cast<size_t>(attr.substr(4)), value);
157 return;
158 }
159 } catch (boost::bad_lexical_cast&) {
160 }
161
162 oself.attr("__class__").attr("__base__").attr("__setattr__")(oself, attr, value);
163 }
164
171 Tensor2<double> activeCond(size_t n, double U, double jy, double T) override {
172 if (n >= this->active.size() || n >= cond_function.size() || cond_function[n].is_none())
173 throw IndexError("no conductivity for active region {}", n);
174 py::object cond = cond_function[n](U, jy, T);
175 py::extract<double> double_cond(cond);
176 if (double_cond.check()) return Tensor2<double>(0., double_cond());
177 return py::extract<Tensor2<double>>(cond);
178 }
179
180 std::string getClassName() const override;
181};
182
183template <typename Class> void setCondJunc(Class& self, const py::object& value) {
184 py::extract<double> double_cond(value);
185 if (double_cond.check())
186 self.setCondJunc(double_cond());
187 else
188 self.setCondJunc(py::extract<Tensor2<double>>(value)());
189}
190
191template <> std::string PythonCondSolver<Geometry2DCartesian>::getClassName() const { return "electrical.ActiveCond2D"; }
192template <> std::string PythonCondSolver<Geometry2DCylindrical>::getClassName() const { return "electrical.ActiveCondCyl"; }
193template <> std::string PythonCondSolver<Geometry3D>::getClassName() const { return "electrical.ActiveCond3D"; }
194
195template <typename __Class__>
196inline static ExportSolver<__Class__> register_electrical_solver(const char* name, const char* geoname) {
197 ExportSolver<__Class__> solver(name,
198 format(u8"{0}(name=\"\")\n\n"
199 u8"Finite element thermal solver for {1} geometry.",
200 name, geoname)
201 .c_str(),
202 py::init<std::string>(py::arg("name") = ""));
203 METHOD(compute, compute, u8"Run electrical calculations", py::arg("loops") = 0);
204 METHOD(get_total_current, getTotalCurrent, u8"Get total current flowing through active region (mA)", py::arg("nact") = 0);
205 RO_PROPERTY(err, getErr, u8"Maximum estimated error");
206 RECEIVER(inTemperature, u8"");
207 PROVIDER(outVoltage, u8"");
208 PROVIDER(outCurrentDensity, u8"");
209 PROVIDER(outHeat, u8"");
210 PROVIDER(outConductivity, u8"");
211 BOUNDARY_CONDITIONS(voltage_boundary, u8"Boundary conditions of the first kind (constant potential)");
212 RW_FIELD(maxerr, u8"Limit for the potential updates");
213 RW_FIELD(convergence, u8"Convergence method.\n\nIf stable, covergence is slown down to ensure stability.");
214 RW_PROPERTY(pcond, getCondPcontact, setCondPcontact, u8"Conductivity of the p-contact");
215 RW_PROPERTY(ncond, getCondNcontact, setCondNcontact, u8"Conductivity of the n-contact");
216 solver.add_property("start_cond", &__Class__::getCondJunc, &setCondJunc<__Class__>,
217 u8"Default effective conductivity of the active region.\n\n"
218 u8"Effective junction conductivity will be computed starting from this value.\n"
219 u8"Note that the actual junction conductivity after convergence can be obtained\n"
220 u8"with :attr:`outConductivity`.");
221 solver.attr("pnjcond") = solver.attr("start_cond");
222 solver.add_property("outPotential", outPotential, u8"Not available in this solver. Use :attr:`outVoltage` instead.");
223 METHOD(get_electrostatic_energy, getTotalEnergy,
224 u8"Get the energy stored in the electrostatic field in the analyzed structure.\n\n"
225 u8"Return:\n"
226 u8" Total electrostatic energy [J].\n");
227 METHOD(get_capacitance, getCapacitance,
228 u8"Get the structure capacitance.\n\n"
229 u8"Return:\n"
230 u8" Total capacitance [pF].\n\n"
231 u8"Note:\n"
232 u8" This method can only be used it there are exactly two boundary conditions\n"
233 u8" specifying the voltage. Otherwise use :meth:`get_electrostatic_energy` to\n"
234 u8" obtain the stored energy $W$ and compute the capacitance as:\n"
235 u8" $C = 2 \\, W / U^2$, where $U$ is the applied voltage.\n");
236 METHOD(get_total_heat, getTotalHeat,
237 u8"Get the total heat produced by the current flowing in the structure.\n\n"
238 u8"Return:\n"
239 u8" Total produced heat (mW).\n");
241
242 return solver;
243}
244
245template <typename GeoT> inline static void register_shockley_solver(const char* name, const char* geoname) {
248 solver.add_property("beta", &__Class__::getBeta0, &__Class__::setBeta0,
249 u8"Junction coefficient (1/V).\n\n"
250 u8"In case, there is more than one junction you may set $\\\\beta$ parameter for any\n"
251 u8"of them by using ``beta#`` property, where # is the junction number (specified\n"
252 u8"by a role ``junction#`` or ``active#``).\n\n"
253 u8"``beta`` is an alias for ``beta0``.\n");
254 solver.add_property("js", &__Class__::getJs0, &__Class__::setJs0,
255 u8"Reverse bias current density (A/m\\ :sup:`2`\\ ).\n\n"
256 u8"In case, there is more than one junction you may set $j_s$ parameter for any\n"
257 u8"of them by using ``js#`` property, where # is the junction number (specified\n"
258 u8"by a role ``junction#`` or ``active#``).\n\n"
259 u8"``js`` is an alias for ``js0``.\n");
260 solver.def("__getattr__", &__Class__::__getattr__);
261 solver.def("__setattr__", &__Class__::__setattr__);
262}
263
264template <typename GeoT> inline static void register_cond_solver(const char* name, const char* geoname) {
267 solver.add_property("cond", &__Class__::getCond0, &__Class__::setCond0,
268 u8"Junction conductivity function [S/m].\n\n"
269 u8"This function should take junction voltage (V), current density (kA/cm²)\n"
270 u8"and temperature (K) as arguments and return a conductivity [S/m]. In case,\n"
271 u8"there is more than one junction you may set such function for any of them\n"
272 u8"by using ``cond#`` property, where # is the junction number (specified by\n"
273 u8"a role ``junction#`` or ``active#``).\n\n"
274 u8"Example:\n\n"
275 u8" >>> def cond(U, j, T):\n"
276 u8" ... return 1e-3 * j**2 + 1e-4 * T**2\n"
277 u8" >>> solver.cond = cond\n\n"
278 u8"``cond`` is an alias for ``cond0``.\n");
279 solver.def("__getattr__", &__Class__::__getattr__);
280 solver.def("__setattr__", &__Class__::__setattr__);
281}
282
290 register_shockley_solver<Geometry2DCartesian>("Shockley2D", "2D Cartesian");
291 register_shockley_solver<Geometry2DCylindrical>("ShockleyCyl", "2D cylindrical");
292 register_shockley_solver<Geometry3D>("Shockley3D", "3D Cartesian");
293
294 register_cond_solver<Geometry2DCartesian>("ActiveCond2D", "2D Cartesian");
295 register_cond_solver<Geometry2DCylindrical>("ActiveCondCyl", "2D cylindrical");
296 register_cond_solver<Geometry3D>("ActiveCond3D", "3D Cartesian");
297}