|
PLaSK library
|
Solver performing calculations in 2D Cartesian or Cylindrical space using finite element method. More...
#include <solvers/electrical/shockley/electr2d.hpp>
Classes | |
| struct | Active |
| Details of active region. More... | |
Public Member Functions | |
| double | compute (unsigned loops=1) |
| Run electrical calculations. | |
| double | integrateCurrent (size_t vindex, bool onlyactive=false) |
| Integrate vertical total current at certain level. | |
| double | getTotalCurrent (size_t nact=0) |
| Integrate vertical total current flowing vertically through active region. | |
| double | getTotalEnergy () |
| Compute total electrostatic energy stored in the structure. | |
| double | getCapacitance () |
| Estimate structure capacitance. | |
| double | getTotalHeat () |
| Compute total heat generated by the structure in unit time. | |
| double | getErr () const |
| Return the maximum estimated error. | |
| double | getCondPcontact () const |
| Get p-contact layer conductivity [S/m]. | |
| void | setCondPcontact (double cond) |
| Set p-contact layer conductivity [S/m]. | |
| double | getCondNcontact () const |
| Get n-contact layer conductivity [S/m]. | |
| void | setCondNcontact (double cond) |
| Set n-contact layer conductivity [S/m]. | |
| DataVector< const Tensor2< double > > | getCondJunc () const |
| Get data with junction effective conductivity. | |
| void | setCondJunc (double cond) |
| Set junction effective conductivity to the single value. | |
| void | setCondJunc (Tensor2< double > cond) |
| Set junction effective conductivity to the single value. | |
| void | setCondJunc (const DataVector< Tensor2< double > > &cond) |
| Set junction effective conductivity to previously read data. | |
| void | loadConfiguration (XMLReader &source, Manager &manager) override |
Load configuration from given source. | |
| void | parseConfiguration (XMLReader &source, Manager &manager) |
| ElectricalFem2DSolver (const std::string &name="") | |
| ~ElectricalFem2DSolver () | |
| double | integrateCurrent (size_t vindex, bool onlyactive) |
| double | integrateCurrent (size_t vindex, bool onlyactive) |
| double | getTotalEnergy () |
| double | getTotalEnergy () |
| double | getTotalHeat () |
| double | getTotalHeat () |
Public Member Functions inherited from plask::FemSolverWithMaskedMesh< Geometry2DType, RectangularMesh< 2 > > | |
| FemSolverWithMaskedMesh (const std::string &name="") | |
| EmptyElementsHandling | getEmptyElements () const |
| Are we using full mesh? | |
| void | setEmptyElements (EmptyElementsHandling val) |
| Set whether we should use full mesh. | |
| bool | parseFemConfiguration (XMLReader &reader, Manager &manager) |
| void | setupMaskedMesh () |
| void | onInitialize () |
| Initialize the solver. | |
| FemMatrix * | getMatrix () |
| FemMatrix * | getMatrix () |
Public Member Functions inherited from plask::FemSolverWithMesh< SpaceT, MeshT > | |
| FemSolverWithMesh (const std::string &name="") | |
| bool | parseFemConfiguration (XMLReader &reader, Manager &manager) |
| FemMatrix * | getMatrix () |
| FemMatrix * | getMatrix () |
Public Member Functions inherited from plask::SolverWithMesh< SpaceT, MeshT > | |
| SolverWithMesh (const std::string &name="") | |
| ~SolverWithMesh () | |
| void | loadConfiguration (XMLReader &source, Manager &manager) override |
Load configuration from given source. | |
| void | parseStandardConfiguration (XMLReader &source, Manager &manager, const std::string &expected_msg="solver configuration element") |
| virtual void | onMeshChange (const typename MeshT::Event &PLASK_UNUSED(evt)) |
| This method is called just after the mesh has been changed. | |
| void | onGeometryChange (const Geometry::Event &PLASK_UNUSED(evt)) override |
| This method is called when the geometry is changed. | |
| MeshT & | meshRef () const |
| Get current module mesh. | |
| shared_ptr< MeshT > | getMesh () const |
| Get current solver mesh. | |
| void | setMesh (const shared_ptr< MeshT > &mesh) |
| Set new mesh for the solver. | |
| void | setMesh (shared_ptr< MeshGeneratorD< MeshT::DIM > > generator) |
| Set new mesh got from generator. | |
Public Member Functions inherited from plask::SolverOver< SpaceT > | |
| SolverOver (const std::string &name="") | |
| ~SolverOver () | |
| void | parseStandardConfiguration (XMLReader &source, Manager &manager, const std::string &expected_msg="solver configuration element") |
| virtual void | onGeometryChange (const Geometry::Event &) |
| This method is called when the geometry is changed. | |
| shared_ptr< SpaceT > | getGeometry () const |
| Get current solver geometry space. | |
| void | setGeometry (const shared_ptr< SpaceT > &geometry) |
| Set new geometry for the solver. | |
Public Member Functions inherited from plask::Solver | |
| bool | initCalculation () |
| This should be called on beginning of each calculation method to ensure that solver will be initialized. | |
| Solver (const std::string &name="") | |
| Construct uninitialized solver. | |
| virtual | ~Solver () |
| Virtual destructor (for subclassing). Do nothing. | |
| void | parseStandardConfiguration (XMLReader &source, Manager &manager, const std::string &expected_msg="solver configuration element") |
Load standard configuration (geometry, mesh) tags from source. | |
| bool | isInitialized () |
| Check if solver is already initialized. | |
| void | invalidate () |
| This method should be and is called if something important was changed: calculation space, mesh, etc. | |
| virtual std::string | getClassName () const =0 |
| Get name of solver. | |
| std::string | getId () const |
| Get solver id. | |
| std::string | getName () const |
| virtual std::string | getClassDescription () const |
| Get a description of this solver. | |
| template<typename ArgT = double, typename ValT = double> | |
| DataLog< ArgT, ValT > | dataLog (const std::string &chart_name, const std::string &axis_arg_name, const std::string &axis_val_name) |
| template<typename ArgT = double, typename ValT = double> | |
| DataLog< ArgT, ValT > | dataLog (const std::string &axis_arg_name, const std::string &axis_val_name) |
| template<typename ... Args> | |
| void | writelog (LogLevel level, std::string msg, Args &&... params) const |
| Log a message for this solver. | |
Public Attributes | |
| double | maxerr |
| Maximum relative current density correction accepted as convergence. | |
| BoundaryConditions< RectangularMesh< 2 >::Boundary, double > | voltage_boundary |
| Boundary condition. | |
| ProviderFor< Voltage, Geometry2DType >::Delegate | outVoltage |
| ProviderFor< CurrentDensity, Geometry2DType >::Delegate | outCurrentDensity |
| ProviderFor< Heat, Geometry2DType >::Delegate | outHeat |
| ProviderFor< Conductivity, Geometry2DType >::Delegate | outConductivity |
| ReceiverFor< Temperature, Geometry2DType > | inTemperature |
| Convergence | convergence |
| Convergence method. | |
Public Attributes inherited from plask::FemSolverWithMesh< SpaceT, MeshT > | |
| FemMatrixAlgorithm | algorithm = ALGORITHM_CHOLESKY |
| Factorization algorithm to use. | |
| IterativeMatrixParams | iter_params |
| Parameters of iterative solver. | |
Protected Member Functions | |
| void | setLocalMatrix (double &k44, double &k33, double &k22, double &k11, double &k43, double &k21, double &k42, double &k31, double &k32, double &k41, double ky, double width, const Vec< 2, double > &midpoint) |
| Save locate stiffness matrix to global one. | |
| virtual Tensor2< double > | activeCond (size_t n, double U, double jy, double T)=0 |
| Compute conductivity int the the active region. | |
| LazyData< double > | loadConductivities () |
| Load conductivities. | |
| void | saveConductivities () |
| Save conductivities of active region. | |
| void | saveHeatDensities () |
| Create 2D-vector with calculated heat densities. | |
| void | onInitialize () override |
| Initialize the solver. | |
| void | onInvalidate () override |
| Invalidate the data. | |
| void | setupActiveRegions () |
| Get info on active region. | |
| void | setMatrix (FemMatrix &A, DataVector< double > &B, const BoundaryConditionsWithMesh< RectangularMesh< 2 >::Boundary, double > &bvoltage, const LazyData< double > &temperature) |
| Set stiffness matrix + load vector. | |
| size_t | isActive (const Vec< 2 > &point) const |
Return true if the specified point is at junction. | |
| size_t | isActive (const RectangularMaskedMesh2D::Element &element) const |
Return true if the specified element is a junction. | |
| const LazyData< double > | getVoltage (shared_ptr< const MeshD< 2 > > dest_mesh, InterpolationMethod method) const |
| const LazyData< double > | getHeatDensities (shared_ptr< const MeshD< 2 > > dest_mesh, InterpolationMethod method) |
| const LazyData< Vec< 2 > > | getCurrentDensities (shared_ptr< const MeshD< 2 > > dest_mesh, InterpolationMethod method) |
| const LazyData< Tensor2< double > > | getConductivity (shared_ptr< const MeshD< 2 > > dest_mesh, InterpolationMethod method) |
| void | setLocalMatrix (double &, double &, double &, double &, double &, double &, double &, double &, double &, double &, double, double, const Vec< 2, double > &) |
| void | setLocalMatrix (double &k44, double &k33, double &k22, double &k11, double &k43, double &k21, double &k42, double &k31, double &k32, double &k41, double, double, const Vec< 2, double > &midpoint) |
Protected Member Functions inherited from plask::SolverOver< SpaceT > | |
| template<typename Boundary , typename ConditionT > | |
| void | readBoundaryConditions (Manager &manager, XMLReader &reader, BoundaryConditions< Boundary, ConditionT > &dest) |
| Read boundary conditions using information about the geometry of this solver. | |
Protected Attributes | |
| double | pcond |
| p-contact electrical conductivity [S/m] | |
| double | ncond |
| n-contact electrical conductivity [S/m] | |
| int | loopno |
| Number of completed loops. | |
| double | toterr |
| Maximum estimated error during all iterations (useful for single calculations managed by external python script) | |
| Vec< 2, double > | maxcur |
| Maximum current in the structure. | |
| DataVector< Tensor2< double > > | junction_conductivity |
| electrical conductivity for p-n junction in y-direction [S/m] | |
| Tensor2< double > | default_junction_conductivity |
| default electrical conductivity for p-n junction in y-direction [S/m] | |
| DataVector< Tensor2< double > > | conds |
| Cached element conductivities. | |
| DataVector< double > | potentials |
| Computed potentials. | |
| DataVector< Vec< 2, double > > | currents |
| Computed current densities. | |
| DataVector< double > | heats |
| Computed and cached heat source densities. | |
| std::vector< Active > | active |
| Active regions information. | |
Protected Attributes inherited from plask::FemSolverWithMaskedMesh< Geometry2DType, RectangularMesh< 2 > > | |
| plask::shared_ptr< RectangularMaskedMesh< MeshT::DIM > > | maskedMesh |
| EmptyElementsHandling | empty_elements |
| Should we use full mesh? | |
Protected Attributes inherited from plask::SolverWithMesh< SpaceT, MeshT > | |
| shared_ptr< MeshT > | mesh |
| Mesh over which the calculations are performed. | |
| boost::signals2::connection | mesh_signal_connection |
| Connection of mesh to onMeshChange method, see http://www.boost.org/doc/libs/1_55_0/doc/html/signals2/tutorial.html#idp204830936. | |
Protected Attributes inherited from plask::SolverOver< SpaceT > | |
| shared_ptr< SpaceT > | geometry |
| Space in which the calculations are performed. | |
Protected Attributes inherited from plask::Solver | |
| bool | initialized |
true only if solver is initialized | |
Additional Inherited Members | |
Public Types inherited from plask::SolverWithMesh< SpaceT, MeshT > | |
| typedef MeshT | MeshType |
| Type of the mesh for this solver. | |
Public Types inherited from plask::SolverOver< SpaceT > | |
| typedef SpaceT | SpaceType |
| of the space for this solver | |
Solver performing calculations in 2D Cartesian or Cylindrical space using finite element method.
Definition at line 25 of file electr2d.hpp.
| plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::ElectricalFem2DSolver | ( | const std::string & | name = "" | ) |
Definition at line 19 of file electr2d.cpp.
| plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::~ElectricalFem2DSolver | ( | ) |
Definition at line 88 of file electr2d.cpp.
|
protectedpure virtual |
Compute conductivity int the the active region.
| n | active region number |
| U | junction voltage (V) |
| jy | vertical current (kA/cm²) |
| T | temperature (K) |
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::compute | ( | unsigned | loops = 1 | ) |
Run electrical calculations.
| loops | maximum number of loops to run |
Definition at line 363 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::getCapacitance | ( | ) |
Estimate structure capacitance.
Definition at line 617 of file electr2d.cpp.
|
inline |
Get data with junction effective conductivity.
Definition at line 238 of file electr2d.hpp.
|
inline |
Get n-contact layer conductivity [S/m].
Definition at line 230 of file electr2d.hpp.
|
inline |
Get p-contact layer conductivity [S/m].
Definition at line 222 of file electr2d.hpp.
|
protected |
Definition at line 565 of file electr2d.cpp.
|
protected |
Definition at line 520 of file electr2d.cpp.
|
inline |
Return the maximum estimated error.
Definition at line 219 of file electr2d.hpp.
|
protected |
Definition at line 542 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::getTotalCurrent | ( | size_t | nact = 0 | ) |
Integrate vertical total current flowing vertically through active region.
| nact | number of the active region |
Definition at line 498 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DCartesian >::getTotalEnergy | ( | ) |
Definition at line 574 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DCylindrical >::getTotalEnergy | ( | ) |
Definition at line 595 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::getTotalEnergy | ( | ) |
Compute total electrostatic energy stored in the structure.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DCartesian >::getTotalHeat | ( | ) |
Definition at line 627 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DCylindrical >::getTotalHeat | ( | ) |
Definition at line 639 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::getTotalHeat | ( | ) |
Compute total heat generated by the structure in unit time.
|
protected |
Definition at line 508 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DCartesian >::integrateCurrent | ( | size_t | vindex, |
| bool | onlyactive | ||
| ) |
Definition at line 466 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DCylindrical >::integrateCurrent | ( | size_t | vindex, |
| bool | onlyactive | ||
| ) |
Definition at line 481 of file electr2d.cpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::integrateCurrent | ( | size_t | vindex, |
| bool | onlyactive = false |
||
| ) |
Integrate vertical total current at certain level.
| vindex | vertical index of the element mesh to perform integration at |
| onlyactive | if true only current in the active region is considered |
|
inlineprotected |
Return true if the specified element is a junction.
Definition at line 158 of file electr2d.hpp.
|
inlineprotected |
Return true if the specified point is at junction.
| point | point to test |
Definition at line 132 of file electr2d.hpp.
|
protected |
|
overridevirtual |
Load configuration from given source.
XML reader (source) point to opening of this solver tag and after return from this method should point to this solver closing tag.
| source | source of configuration |
| manager | manager from which information about geometry, meshes, materials, and so on can be get if needed |
Reimplemented from plask::Solver.
Definition at line 37 of file electr2d.cpp.
|
overrideprotectedvirtual |
Initialize the solver.
Reimplemented from plask::Solver.
Definition at line 171 of file electr2d.cpp.
|
overrideprotectedvirtual |
| void plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::parseConfiguration | ( | XMLReader & | source, |
| Manager & | manager | ||
| ) |
Definition at line 42 of file electr2d.cpp.
|
protected |
Save conductivities of active region.
Definition at line 355 of file electr2d.cpp.
|
protected |
Create 2D-vector with calculated heat densities.
Definition at line 442 of file electr2d.cpp.
|
inline |
Set junction effective conductivity to previously read data.
Definition at line 250 of file electr2d.hpp.
|
inline |
Set junction effective conductivity to the single value.
Definition at line 240 of file electr2d.hpp.
|
inline |
Set junction effective conductivity to the single value.
Definition at line 245 of file electr2d.hpp.
|
inline |
Set n-contact layer conductivity [S/m].
Definition at line 232 of file electr2d.hpp.
|
inline |
Set p-contact layer conductivity [S/m].
Definition at line 224 of file electr2d.hpp.
|
inlineprotected |
Definition at line 190 of file electr2d.cpp.
|
inlineprotected |
Save locate stiffness matrix to global one.
|
inlineprotected |
Definition at line 207 of file electr2d.cpp.
|
protected |
Set stiffness matrix + load vector.
Definition at line 235 of file electr2d.cpp.
|
protected |
Get info on active region.
Definition at line 90 of file electr2d.cpp.
|
protected |
Active regions information.
Definition at line 67 of file electr2d.hpp.
|
protected |
Cached element conductivities.
Definition at line 62 of file electr2d.hpp.
| Convergence plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::convergence |
Convergence method.
Definition at line 176 of file electr2d.hpp.
|
protected |
Computed current densities.
Definition at line 64 of file electr2d.hpp.
|
protected |
default electrical conductivity for p-n junction in y-direction [S/m]
Definition at line 60 of file electr2d.hpp.
|
protected |
Computed and cached heat source densities.
Definition at line 65 of file electr2d.hpp.
| ReceiverFor<Temperature, Geometry2DType> plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::inTemperature |
Definition at line 174 of file electr2d.hpp.
|
protected |
electrical conductivity for p-n junction in y-direction [S/m]
Definition at line 59 of file electr2d.hpp.
|
protected |
Number of completed loops.
Definition at line 54 of file electr2d.hpp.
|
protected |
Maximum current in the structure.
Definition at line 57 of file electr2d.hpp.
| double plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::maxerr |
Maximum relative current density correction accepted as convergence.
Definition at line 161 of file electr2d.hpp.
|
protected |
n-contact electrical conductivity [S/m]
Definition at line 52 of file electr2d.hpp.
| ProviderFor<Conductivity,Geometry2DType>::Delegate plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::outConductivity |
Definition at line 172 of file electr2d.hpp.
| ProviderFor<CurrentDensity,Geometry2DType>::Delegate plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::outCurrentDensity |
Definition at line 168 of file electr2d.hpp.
| ProviderFor<Heat,Geometry2DType>::Delegate plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::outHeat |
Definition at line 170 of file electr2d.hpp.
| ProviderFor<Voltage,Geometry2DType>::Delegate plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::outVoltage |
Definition at line 166 of file electr2d.hpp.
|
protected |
p-contact electrical conductivity [S/m]
Definition at line 51 of file electr2d.hpp.
|
protected |
Computed potentials.
Definition at line 63 of file electr2d.hpp.
|
protected |
Maximum estimated error during all iterations (useful for single calculations managed by external python script)
Definition at line 55 of file electr2d.hpp.
| BoundaryConditions<RectangularMesh<2>::Boundary, double> plask::electrical::shockley::ElectricalFem2DSolver< Geometry2DType >::voltage_boundary |
Boundary condition.
Definition at line 164 of file electr2d.hpp.