PLaSK library
Loading...
Searching...
No Matches
thermoelectric.py
Go to the documentation of this file.
1# This file is part of PLaSK (https://plask.app) by Photonics Group at TUL
2# Copyright (c) 2022 Lodz University of Technology
3#
4# This program is free software: you can redistribute it and/or modify
5# it under the terms of the GNU General Public License as published by
6# the Free Software Foundation, version 3.
7#
8# This program is distributed in the hope that it will be useful,
9# but WITHOUT ANY WARRANTY; without even the implied warranty of
10# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11# GNU General Public License for more details.
12
13# coding: utf8
14# Copyright (C) 2014 Photonics Group, Lodz University of Technology
15
16import plask
17
18import thermal.static
19import electrical.shockley
20
21
22class attribute(object):
23 def __init__(self, value):
24 self.value = value
25 def __repr__(self):
26 return str(self.value)
27
28
29# Get unique suffix for savefiles
30if plask.BATCH:
31 suffix = "-{}".format(plask.JOBID)
32else:
33 import time as _time
34 suffix = _time.strftime("-%Y%m%d-%H%M", _time.localtime(plask.JOBID))
35
36
37def h5open(filename, group):
38 import h5py
39 if filename is None:
40 import sys
41 from os.path import basename, abspath
42 filename = sys.argv[0]
43 if filename.endswith('.py'): filename = filename[:-3]
44 elif filename.endswith('.xpl'): filename = filename[:-4]
45 elif filename == '': filename = 'console'
46 filename = abspath(basename(filename))
47 filename += suffix + '.h5'
48 if type(filename) is h5py.File:
49 h5file = filename
50 filename = h5file.filename
51 close = False
52 else:
53 h5file = h5py.File(filename, 'a')
54 close = True
55 group_parts = [g for g in group.split('/') if g]
56 orig_group = group = '/'.join(group_parts)
57 idx = 1
58 while group in h5file:
59 group = '/'.join(["{}-{}".format(group_parts[0], idx)] + group_parts[1:])
60 idx += 1
61 if group is not orig_group:
62 plask.print_log('warning',
63 "Group '/{}' exists in HDF5 file '{}'. Saving to group '/{}'" \
64 .format(orig_group, filename, group))
65 group = '/' + group
66 return h5file, group, filename, close
67
68
70
71 _Thermal = None
72 _Electrical = None
73
74 def __init__(self, name):
75 super().__init__(name)
76 self.thermal = self._Thermal(name)
77 self.electrical = self._Electrical(name)
78 self.__reconnect()
79
80 self.tfreq = 6
81
82 def __reconnect(self):
83 self.electrical.inTemperature = self.thermal.outTemperature
84 self.thermal.inHeat = self.electrical.outHeat
85
86 def reconnect(self):
87 """
88 Reconnect all internal solvers.
89
90 This method should be called if some of the internal solvers were changed manually.
91 """
92 self.__reconnect()
93
94 def _read_attr(self, tag, attr, solver, type=None, pyattr=None):
95 if pyattr is None: pyattr = attr
96 if attr in tag:
97 try:
98 if type is not None:
99 val = type(tag[attr])
100 else:
101 val = tag[attr]
102 except ValueError:
103 raise plask.XMLError("{}: {} attribute {} has illegal value '{}'".format(
104 tag, type.__name__.title(), attr, tag[attr]))
105 else:
106 setattr(solver, pyattr, val)
107
108 def load_xpl(self, xpl, manager):
109 for tag in xpl:
110 self._parse_xpl(tag, manager)
111
112 def _parse_fem_tag(self, tag, manager, solver):
113 self._read_attr(tag, 'algorithm', solver)
114 for it in tag:
115 if it == 'iterative':
116 self._read_attr(it, 'accelerator', solver.iterative)
117 self._read_attr(it, 'preconditioner', solver.iterative)
118 self._read_attr(it, 'noconv', solver.iterative)
119 self._read_attr(it, 'maxit', solver.iterative, int)
120 self._read_attr(it, 'maxerr', solver.iterative, float)
121 self._read_attr(it, 'nfact', solver.iterative, int)
122 self._read_attr(it, 'omega', solver.iterative, float)
123 self._read_attr(it, 'ndeg', solver.iterative, int)
124 self._read_attr(it, 'lvfill', solver.iterative, int)
125 self._read_attr(it, 'ltrunc', solver.iterative, int)
126 self._read_attr(it, 'nsave', solver.iterative, int)
127 self._read_attr(it, 'nrestart', solver.iterative, int)
128 else:
129 raise plask.XMLError("{}: Unrecognized tag '{}'".format(it, it.name))
130
131 def _parse_xpl(self, tag, manager):
132 if tag == 'geometry':
133 self.thermal.geometry = tag.getitem(manager.geo, 'thermal')
134 self.electrical.geometry = tag.getitem(manager.geo, 'electrical')
135 elif tag == 'mesh':
136 self.thermal.mesh = tag.getitem(manager.msh, 'thermal')
137 self._read_attr(tag, 'empty-elements', self.electrical, str, 'empty_elements')
138 self.electrical.mesh = tag.getitem(manager.msh, 'electrical')
139 elif tag == 'junction':
140 self._read_attr(tag, 'pnjcond', self.electrical, float)
141 for key, val in tag.attrs.items():
142 if key.startswith('beta') or key.startswith('js'):
143 try:
144 val = float(val)
145 except ValueError:
146 raise plask.XMLError("{}: Float attribute '{}' has illegal value '{}'".format(tag, key, val))
147 setattr(self.electrical, key, val)
148 elif tag == 'contacts':
149 self._read_attr(tag, 'pcond', self.electrical, float)
150 self._read_attr(tag, 'ncond', self.electrical, float)
151 elif tag == 'loop':
152 self.tfreq = int(tag.get('tfreq', self.tfreq))
153 self._read_attr(tag, 'inittemp', self.thermal, float, 'inittemp')
154 self._read_attr(tag, 'maxterr', self.thermal, float, 'maxerr')
155 self._read_attr(tag, 'maxcerr', self.electrical, float, 'maxerr')
156 elif tag == 'tmatrix':
157 self._parse_fem_tag(tag, manager, self.thermal)
158 elif tag == 'ematrix':
159 self._parse_fem_tag(tag, manager, self.electrical)
160 elif tag == 'temperature':
161 self.thermal.temperature_boundary.read_from_xpl(tag, manager)
162 elif tag == 'heatflux':
163 self.thermal.heatflux_boundary.read_from_xpl(tag, manager)
164 elif tag == 'convection':
165 self.thermal.convection_boundary.read_from_xpl(tag, manager)
166 elif tag == 'radiation':
167 self.thermal.radiation_boundary.read_from_xpl(tag, manager)
168 elif tag == 'voltage':
169 self.electrical.voltage_boundary.read_from_xpl(tag, manager)
170 else:
171 raise plask.XMLError("{}: Unrecognized tag '{}'".format(tag, tag.name))
172
173 def on_initialize(self):
174 self.thermal.initialize()
175 self.electrical.initialize()
176
177 def on_invalidate(self):
178 self.thermal.invalidate()
180
181 def compute(self, save=True, invalidate=True, group='ThermoElectric'):
182 """
183 Run calculations.
184
185 In the beginning the solvers are invalidated and next, the thermo-
186 electric algorithm is executed until both solvers converge to the
187 value specified in their configuration in the `maxerr` property.
188
189 Args:
190 save (bool or str): If `True` the computed fields are saved to the
191 HDF5 file named after the script name with the suffix denoting
192 either the batch job id or the current time if no batch system
193 is used. The filename can be overridden by setting this parameter
194 as a string.
195
196 invalidate (bool): If this flag is set, solvers are invalidated
197 in the beginning of the computations.
198
199 group (str): HDF5 group to save the data under.
200 """
201 if invalidate:
202 self.thermal.invalidate()
204
205 self.initialize()
206
207 verr = 2. * self.electrical.maxerr
208 terr = 2. * self.thermal.maxerr
209 while terr > self.thermal.maxerr or verr > self.electrical.maxerr:
210 verr = self.electrical.compute(self.tfreq)
211 terr = self.thermal.compute(1)
212
213 infolines = self._get_info()
214 plask.print_log('important', "Thermoelectric Computations Finished")
215 for line in infolines:
216 plask.print_log('important', " " + line)
217
218 if save:
219 self.save(None if save is True else save, group)
220
221 def get_total_current(self, nact=0):
222 """
223 Get total current flowing through active region (mA)
224 """
225 return self.electrical.get_total_current(nact)
226
227 @staticmethod
228 def _iter_levels(geometry, mesh, *required):
229 if isinstance(mesh, (plask.mesh.Mesh1D, plask.ndarray, list, tuple)):
230 hor = mesh,
231 Mesh = plask.mesh.Rectangular2D
232 elif isinstance(mesh, plask.mesh.Rectangular2D):
233 hor = mesh.axis0,
234 Mesh = plask.mesh.Rectangular2D
235 elif isinstance(mesh, plask.mesh.Rectangular3D):
236 hor = mesh.axis0, mesh.axis1
237 Mesh = plask.mesh.Rectangular3D
238 else:
239 return
240
241 points = Mesh.SimpleGenerator()(geometry).elements.mesh
242 levels = {}
243 for p in points:
244 roles = geometry.get_roles(p)
245 suffix = None
246 for role in roles:
247 for prefix in ('active', 'junction'):
248 if role.startswith(prefix):
249 suffix = role[len(prefix):]
250 break
251 if suffix is not None and (not required or any(role in required for role in roles)):
252 levels[suffix] = p[-1]
253
254 for name, v in levels.items():
255 axs = hor + ([v],)
256 mesh2 = Mesh(*axs)
257 yield name, mesh2
258
259 def _save_thermoelectric(self, h5file, group):
260 tmesh = self.thermal.mesh
261 vmesh = self.electrical.mesh
262 jmesh = vmesh.elements.mesh
263 temp = self.thermal.outTemperature(tmesh)
264 volt = self.electrical.outVoltage(vmesh)
265 curr = self.electrical.outCurrentDensity(jmesh)
266 plask.save_field(temp, h5file, group + '/Temperature')
267 plask.save_field(volt, h5file, group + '/Potential')
268 plask.save_field(curr, h5file, group + '/CurrentDensity')
269 for name, jmesh2 in self._iter_levels(self.electrical.geometry, jmesh):
270 curr2 = self.electrical.outCurrentDensity(jmesh2)
271 plask.save_field(curr2, h5file, group + '/Junction' + name + 'CurrentDensity')
272
273 def save(self, filename=None, group='ThermoElectric'):
274 """
275 Save the computation results to the HDF5 file.
276
277 Args:
278 filename (str): The file name to save to.
279 If omitted, the file name is generated automatically based on
280 the script name with suffix denoting either the batch job id or
281 the current time if no batch system is used.
282
283 group (str): HDF5 group to save the data under.
284 """
285 h5file, group, filename, close = h5open(filename, group)
286 self._save_thermoelectric(h5file, group)
287 if close:
288 h5file.close()
289 plask.print_log('info', "Fields saved to file '{}' in group '{}'".format(filename, group))
290 return filename
291
293 """
294 Get temperature on a thermal mesh.
295 """
296 return self.thermal.outTemperature(self.thermal.mesh)
297
298 def get_voltage(self):
299 """
300 Get voltage on an electrical mesh.
301 """
302 return self.electrical.outVoltage(self.electrical.mesh)
303
304 def get_vertical_voltage(self, at=0):
305 """
306 Get computed voltage along the vertical axis.
307
308 Args:
309 at (float): Horizontal position of the axis at which the voltage
310 is plotted.
311 """
312 if isinstance(self.electrical.geometry, plask.geometry.Cartesian3D):
313 try:
314 at0, at1 = at
315 except TypeError:
316 at0 = at1 = at
317 mesh = plask.mesh.Rectangular2D(plask.mesh.Ordered([at0]), plask.mesh.Ordered([at1]),
318 self.electrical.mesh.axis2)
319 else:
320 mesh = plask.mesh.Rectangular2D(plask.mesh.Ordered([at]), self.electrical.mesh.axis1)
321 field = self.electrical.outVoltage(mesh)
322 return field
323
324 def get_junction_currents(self, refine=16, interpolation='linear'):
325 """
326 Get current densities at the active regions.
327
328 Args:
329 refine (int): Number of points in the plot between each two points
330 in the computational mesh.
331
332 interpolation (str): Interpolation used when retrieving current density.
333
334 Return:
335 dict: Dictionary of junction current density data.
336 Keys are the junction number.
337 """
338 axis = plask.concatenate([
339 plask.linspace(x, self.electrical.mesh.axis0[i+1], refine+1)
340 for i,x in enumerate(list(self.electrical.mesh.axis0)[:-1])
341 ])
342
343 result = {}
344 for lb, msh in self._iter_levels(self.electrical.geometry, axis):
345 if not lb:
346 lb = 0
347 else:
348 try: lb = int(lb)
349 except ValueError: pass
350 result[lb] = self.electrical.outCurrentDensity(msh, interpolation).array[:,0,1]
351 return result
352
353 def plot_temperature(self, geometry_color='0.75', mesh_color=None, geometry_alpha=0.35, mesh_alpha=0.15,
354 geometry_lw=1.0, mesh_lw=1.0, **kwargs):
355 """
356 Plot computed temperature to the current axes.
357
358 Args:
359 geometry_color (str or ``None``): Matplotlib color specification for
360 the geometry. If ``None``, structure is not plotted.
361
362 mesh_color (str or ``None``): Matplotlib color specification for
363 the mesh. If ``None``, the mesh is not plotted.
364
365 geometry_alpha (float): Geometry opacity (1 — fully opaque, 0 – invisible).
366
367 mesh_alpha (float): Mesh opacity (1 — fully opaque, 0 – invisible).
368
369 geometry_lw (float): Line width for geometry.
370
371 mesh_lw (float): Line width for mesh.
372
373 **kwargs: Keyword arguments passed to the plot function.
374
375 See also:
376 :func:`plask.plot_field` : Plot any field obtained from receivers
377 """
378 field = self.get_temperature()
379 plask.plot_field(field, **kwargs)
380 cbar = plask.colorbar(use_gridspec=True)
381 cbar.set_label("Temperature (K)")
382 if geometry_color is not None:
383 plask.plot_geometry(self.thermal.geometry, color=geometry_color, alpha=geometry_alpha, lw=geometry_lw)
384 if mesh_color is not None:
385 plask.plot_mesh(self.thermal.mesh, color=mesh_color, alpha=mesh_alpha, lw=mesh_lw)
386 plask.window_title("Temperature")
387
388 def plot_voltage(self, geometry_color='0.75', mesh_color=None, geometry_alpha=0.35, mesh_alpha=0.15,
389 geometry_lw=1.0, mesh_lw=1.0, **kwargs):
390 """
391 Plot computed voltage to the current axes.
392
393 Args:
394 geometry_color (str or ``None``): Matplotlib color specification
395 for the geometry. If ``None``, structure is not plotted.
396
397 mesh_color (str or ``None``): Matplotlib color specification for
398 the mesh. If ``None``, the mesh is not plotted.
399
400 geometry_alpha (float): Geometry opacity (1 — fully opaque, 0 – invisible).
401
402 mesh_alpha (float): Mesh opacity (1 — fully opaque, 0 – invisible).
403
404 geometry_lw (float): Line width for geometry.
405
406 mesh_lw (float): Line width for mesh.
407
408 **kwargs: Keyword arguments passed to the :func:`plask.plot_field`.
409
410 See also:
411 :func:`plask.plot_field` : Plot any field obtained from receivers
412 """
413 field = self.get_voltage()
414 plask.plot_field(field, **kwargs)
415 cbar = plask.colorbar(use_gridspec=True)
416 cbar.set_label("Voltage (V)")
417 if geometry_color is not None:
418 plask.plot_geometry(self.electrical.geometry, color=geometry_color, alpha=geometry_alpha, lw=geometry_lw)
419 if mesh_color is not None:
420 plask.plot_mesh(self.electrical.mesh, color=mesh_color, alpha=mesh_alpha, lw=mesh_lw)
421 plask.window_title("Voltage")
422
423 def plot_vertical_voltage(self, at=0., **kwargs):
424 """
425 Plot computed voltage along the vertical axis.
426
427 Args:
428 at (float): Horizontal position of the axis at which the voltage
429 is plotted.
430
431 **kwargs: Keyword arguments passed to the plot function.
432 """
433 field = self.get_vertical_voltage(at)
434 plask.plot(field.mesh.axis1, field, **kwargs)
435 plask.xlabel(u"${}$ (µm)".format(plask.config.axes[-1]))
436 plask.ylabel("Voltage (V)")
437 plask.window_title("Voltage")
438
439 def _plot_hbounds(self, solver):
440 simplemesh = plask.mesh.Rectangular2D.SimpleGenerator()(solver.geometry.item)
441 for x in simplemesh.axis0:
442 plask.axvline(x,
443 linestyle=plask.rc.grid.linestyle, linewidth=plask.rc.grid.linewidth,
444 color=plask.rc.grid.color, alpha=plask.rc.grid.alpha)
445
446 def plot_junction_current(self, refine=16, bounds=True, interpolation='linear', label=None, **kwargs):
447 """
448 Plot current density at the active region.
449
450 Args:
451 refine (int): Number of points in the plot between each two points
452 in the computational mesh.
453
454 bounds (bool): If *True* then the geometry objects boundaries are
455 plotted.
456
457 interpolation (str): Interpolation used when retrieving current density.
458
459 label (str or sequence): Label for each junction. It can be a sequence of
460 consecutive labels for each junction, or a string
461 in which case the same label is used for each
462 junction. If omitted automatic label is generated.
463
464 **kwargs: Keyword arguments passed to the plot function.
465 """
466 axis = plask.concatenate([
467 plask.linspace(x, self.electrical.mesh.axis0[i+1], refine+1)
468 for i,x in enumerate(list(self.electrical.mesh.axis0)[:-1])
469 ])
470
471 i = 0
472 for i, (lb, msh) in enumerate(self._iter_levels(self.electrical.geometry, axis)):
473 curr = self.electrical.outCurrentDensity(msh, interpolation).array[:,0,1]
474 s = sum(curr)
475 if label is None:
476 lab = "Junction {:s}".format(lb)
477 elif isinstance(label, tuple) or isinstance(label, tuple):
478 lab = label[i]
479 else:
480 lab = label
481 plask.plot(msh.axis0, curr if s > 0 else -curr,
482 label=lab, **kwargs)
483 if i > 0:
484 plask.legend(loc='best')
485 plask.xlabel(u"${}$ (µm)".format(plask.config.axes[-2]))
486 plask.ylabel(u"Current Density (kA/cm\xb2)")
487 if bounds:
488 self._plot_hbounds(self.electrical)
489 plask.window_title("Current Density")
490
492 try:
493 import __main__
494 defines = [" {} = {}".format(key, __main__.DEF[key]) for key in __main__.__overrites__]
495 except (NameError, KeyError, AttributeError):
496 defines = []
497 if defines:
498 defines = ["Temporary defines:"] + defines
499 return defines
500
501 def _get_info(self):
502 info = self._get_defines_info()
503 try:
504 info.append("Total current (mA): {:8.3f}".format(self.get_total_current()))
505 except:
506 pass
507 try:
508 info.append("Maximum temperature (K): {:8.3f}".format(max(self.get_temperature())))
509 except:
510 pass
511 return info
512
513
515 """
516 Thermo-electric calculations solver without the optical part.
517
518 This solver performs under-threshold thermo-electrical computations.
519 It computes electric current flow and temperature distribution in a self-
520 consistent loop until desired convergence is reached.
521
522 The computations can be executed using `compute` method, after which
523 the results may be save to the HDF5 file with `save` or presented visually
524 using ``plot_...`` methods. If ``save`` parameter of the :meth:`compute` method
525 is *True* the fields are saved automatically after the computations.
526 The file name is based on the name of the executed script with suffix denoting
527 either the launch time or the identifier of a batch job if a batch system
528 (like SLURM, OpenPBS, or SGE) is used.
529 """
530
531 _Thermal = thermal.static.Static2D
532 _Electrical = electrical.shockley.Shockley2D
533
534 outTemperature = property(lambda self: self.thermalthermal.outTemperature,
535 doc=_Thermal.outTemperature.__doc__)
536
537 outHeatFlux = property(lambda self: self.thermalthermal.outHeatFlux,
538 doc=_Thermal.outHeatFlux.__doc__)
539
540 outThermalConductivity = property(lambda self: self.thermalthermal.outThermalConductivity,
541 doc=_Thermal.outThermalConductivity.__doc__)
542
543 outVoltage = property(lambda self: self.electricalelectrical.outVoltage,
544 doc=_Electrical.outVoltage.__doc__)
545
546 outCurrentDensity = property(lambda self: self.electricalelectrical.outCurrentDensity,
547 doc=_Electrical.outCurrentDensity.__doc__)
548
549 outHeat = property(lambda self: self.electricalelectrical.outHeat,
550 doc=_Electrical.outHeat.__doc__)
551
552 outConductivity = property(lambda self: self.electricalelectrical.outConductivity,
553 doc=_Electrical.outConductivity.__doc__)
554
555 thermal = attribute(_Thermal.__name__+"()")
556 """
557 :class:`thermal.static.Static2D` solver used for thermal calculations.
558 """
559
560 electrical = attribute(_Electrical.__name__+"()")
561 """
562 :class:`electrical.shockley.Shockley2D` solver used for electrical calculations.
563 """
564
565 tfreq = 6.0
566 """
567 Number of electrical iterations per single thermal step.
568
569 As temperature tends to converge faster, it is reasonable to repeat thermal
570 solution less frequently.
571 """
572
573
575 """
576 Thermo-electric calculations solver without the optical part.
577
578 This solver performs under-threshold thermo-electrical computations.
579 It computes electric current flow and temperature distribution in a self-
580 consistent loop until desired convergence is reached.
581
582 The computations can be executed using `compute` method, after which
583 the results may be save to the HDF5 file with `save` or presented visually
584 using ``plot_...`` methods. If ``save`` parameter of the :meth:`compute` method
585 is *True* the fields are saved automatically after the computations.
586 The file name is based on the name of the executed script with suffix denoting
587 either the launch time or the identifier of a batch job if a batch system
588 (like SLURM, OpenPBS, or SGE) is used.
589 """
590
591 _Thermal = thermal.static.StaticCyl
592 _Electrical = electrical.shockley.ShockleyCyl
593
594 outTemperature = property(lambda self: self.thermalthermal.outTemperature,
595 doc=_Thermal.outTemperature.__doc__)
596
597 outHeatFlux = property(lambda self: self.thermalthermal.outHeatFlux,
598 doc=_Thermal.outHeatFlux.__doc__)
599
600 outThermalConductivity = property(lambda self: self.thermalthermal.outThermalConductivity,
601 doc=_Thermal.outThermalConductivity.__doc__)
602
603 outVoltage = property(lambda self: self.electricalelectrical.outVoltage,
604 doc=_Electrical.outVoltage.__doc__)
605
606 outCurrentDensity = property(lambda self: self.electricalelectrical.outCurrentDensity,
607 doc=_Electrical.outCurrentDensity.__doc__)
608
609 outHeat = property(lambda self: self.electricalelectrical.outHeat,
610 doc=_Electrical.outHeat.__doc__)
611
612 outConductivity = property(lambda self: self.electricalelectrical.outConductivity,
613 doc=_Electrical.outConductivity.__doc__)
614
615 thermal = attribute(_Thermal.__name__+"()")
616 """
617 :class:`thermal.static.Static2D` solver used for thermal calculations.
618 """
619
620 electrical = attribute(_Electrical.__name__+"()")
621 """
622 :class:`electrical.shockley.Shockley2D` solver used for electrical calculations.
623 """
624
625 tfreq = 6.0
626 """
627 Number of electrical iterations per single thermal step.
628
629 As temperature tends to converge faster, it is reasonable to repeat thermal
630 solution less frequently.
631 """
632
633
635 """
636 Thermo-electric calculations solver without the optical part.
637
638 This solver performs under-threshold thermo-electrical computations.
639 It computes electric current flow and temperature distribution in a self-
640 consistent loop until desired convergence is reached.
641
642 The computations can be executed using `compute` method, after which
643 the results may be save to the HDF5 file with `save` or presented visually
644 using ``plot_...`` methods. If ``save`` parameter of the :meth:`compute` method
645 is *True* the fields are saved automatically after the computations.
646 The file name is based on the name of the executed script with suffix denoting
647 either the launch time or the identifier of a batch job if a batch system
648 (like SLURM, OpenPBS, or SGE) is used.
649 """
650
651 _Thermal = thermal.static.Static3D
652 _Electrical = electrical.shockley.Shockley3D
653
654 outTemperature = property(lambda self: self.thermalthermal.outTemperature,
655 doc=_Thermal.outTemperature.__doc__)
656
657 outHeatFlux = property(lambda self: self.thermalthermal.outHeatFlux,
658 doc=_Thermal.outHeatFlux.__doc__)
659
660 outThermalConductivity = property(lambda self: self.thermalthermal.outThermalConductivity,
661 doc=_Thermal.outThermalConductivity.__doc__)
662
663 outVoltage = property(lambda self: self.electricalelectrical.outVoltage,
664 doc=_Electrical.outVoltage.__doc__)
665
666 outCurrentDensity = property(lambda self: self.electricalelectrical.outCurrentDensity,
667 doc=_Electrical.outCurrentDensity.__doc__)
668
669 outHeat = property(lambda self: self.electricalelectrical.outHeat,
670 doc=_Electrical.outHeat.__doc__)
671
672 outConductivity = property(lambda self: self.electricalelectrical.outConductivity,
673 doc=_Electrical.outConductivity.__doc__)
674
675 thermal = attribute(_Thermal.__name__+"()")
676 """
677 :class:`thermal.static.Static3D` solver used for thermal calculations.
678 """
679
680 electrical = attribute(_Electrical.__name__+"()")
681 """
682 :class:`electrical.shockley.Shockley3D` solver used for electrical calculations.
683 """
684
685 tfreq = 6.0
686 """
687 Number of electrical iterations per single thermal step.
688
689 As temperature tends to converge faster, it is reasonable to repeat thermal
690 solution less frequently.
691 """
692
693
694__all__ = 'ThermoElectric2D', 'ThermoElectricCyl', 'ThermoElectric3D'