hp2D::Maxwell2D_H Class Referenceabstract
Class for calculating Eddy current problem with Maxwell modell in h formulation. More...
#include <Maxwell2D_H.hh>
Public Types | |
enum | boundaryType { PMC = 0, PEC = 1, MAX_TYPE } |
Boundary type. More... | |
enum | solverType { SUPERLU = 0, SUPERLU2 = 1, BICGSTAB = 2, BICGSTAB2 = 3, BICGSTABSUPERLU = 4 } |
Type of the solver. More... | |
typedef Cmplx | type |
Public Member Functions | |
boundaryType & | bType () |
const boundaryType | bType () const |
Returns boundary type. More... | |
const std::string | bTypeStr () const |
Returns name of boundary type as string. More... | |
virtual Real | dissipation () |
Return dissipation power loss. More... | |
concepts::ElementFormula< concepts::Cmplx2d > * | eField () |
Returns a pointer to the e-Field. More... | |
concepts::ElementFormula< Cmplx > * | hField () |
Returns a pointer to the h-Field (solution + h0). More... | |
virtual Real | magnEnergy () |
Return magnetic energy. More... | |
Maxwell2D_H (concepts::EddyGeometry2D &geom, enum concepts::MaxwellBoundary::boundaryType bType=PMC, enum solverType type=SUPERLU, bool diagPrecond=true, bool afterIter=false, const Real eps=EPS0, const Real omega=OMEGA50, const Real mu=MU0, const uint geomRefAttrib=100) | |
Constructor. More... | |
Maxwell2D_H (concepts::EddyGeometry2D &geom, InputMaxwell2D_H &input, const uint geomRefAttrib=100) | |
void | rebuildMesh (concepts::InputAdaptiveModels &input) |
void | rebuildMesh (const uint l=0, const uint p=1, const uint g=0, const uint subdiv=X|Y) |
Rebuilds only the mesh and sets the polynomial degrees. More... | |
const Vector< Cmplx > * | solution () |
Returns solution vector. More... | |
virtual hpAdaptiveSpaceH1 & | space () const |
Returns the space. More... | |
virtual Space< Real > & | space () const=0 |
Returns the space. More... | |
virtual | ~Maxwell2D_H () |
Protected Types | |
enum | subdivTypes |
Protected Member Functions | |
concepts::SparseMatrix< Real > * | identityMatrix_ (concepts::Space< Real > &spc, concepts::SparseMatrix< Cmplx > *S) |
Calculate identity matrix and add's it to system matrix S . More... | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. More... | |
concepts::SparseMatrix< Cmplx > * | laplaceMatrix_ (concepts::Space< Real > &spc, concepts::SparseMatrix< Cmplx > *S) |
Calculate stiffness matrix and add's it to system matrix S . More... | |
void | linearform_ () |
Calculate the load vector, assumes sigma = 0 inside the coil. More... | |
virtual concepts::Mesh & | mesh_ () |
Mesh. More... | |
virtual const std::string | mshAbbr_ () |
Mesh abbreviation string. More... | |
virtual hpFull & | prebuild_ () |
Space Prebuilder. More... | |
virtual void | rebuildMesh_ (const uint l=0, const uint p=1, const uint g=0, const uint subdiv=X|Y) |
Rebuilds only the mesh and sets the polynomial degrees. More... | |
Protected Attributes | |
bool | afterIter_ |
Nachiteration. More... | |
std::unique_ptr< concepts::BoundaryConditions > | bc_ |
Boundary conditions. More... | |
bool | diagPrecond_ |
Using diagonal preconditioning. More... | |
std::unique_ptr< Real > | dissipation_ |
Dissipation power loss. More... | |
Real | eps_ |
Dielectricity constant. More... | |
concepts::EddyGeometry2D & | geom_ |
Mesh and material constants (sigma, j0) More... | |
const uint | geomRefAttr_ |
Attribute of vertices or edges for geometric refinement. More... | |
concepts::PiecewiseFormulaFun< Cmplx, Real > | iOmegaEps_plus_Sigma_Inv_ |
Piecewise constant formula. More... | |
uint | iterations_ |
Number of iterations for iterative solver. More... | |
std::unique_ptr< Real > | magnEnergy_ |
Magnetic energy. More... | |
double | matrixtime_ |
const Real | mu_ |
Permeability constant. More... | |
Real | omega_ |
Angular frequency. More... | |
const std::string | problemName_ |
Name of the problem. More... | |
std::unique_ptr< concepts::Vector< Cmplx > > | residual_ |
std::unique_ptr< Real > | residualNorm_ |
Euclidian norm of the residual of solving the linear system. More... | |
std::unique_ptr< concepts::Vector< Cmplx > > | rhs_ |
std::unique_ptr< Vector< Cmplx > > | sol_ |
Solution vector. More... | |
double | solvetime_ |
Time to solve the system, to build the matrices, to rebuild the space. More... | |
double | spacetime_ |
bool | statusAfterIter_ |
enum solverType | type_ |
Solver type. More... | |
Private Member Functions | |
void | constructor_ () |
Private constructor. More... | |
void | constructSpace_ () |
Constructs the space. More... | |
void | matrices_ () |
Building the matrices. More... | |
virtual void | solve_ () |
Method for solving, throws exception when it wasn't successfull. More... | |
Private Attributes | |
std::unique_ptr< concepts::SparseMatrix< Cmplx > > | A_ |
Stiffness and System matrix. More... | |
enum boundaryType | bType_ |
Type of boundary condition. More... | |
std::unique_ptr< concepts::SparseMatrix< Real > > | M_ |
Mass matrix. More... | |
std::unique_ptr< concepts::SparseMatrix< Cmplx > > | S_ |
std::unique_ptr< hpAdaptiveSpaceH1 > | spc_ |
Space. More... | |
Friends | |
class | concepts::ModelControl< Maxwell2D_H > |
Detailed Description
Class for calculating Eddy current problem with Maxwell modell in h formulation.
Definition at line 122 of file Maxwell2D_H.hh.
Member Typedef Documentation
◆ type
|
inherited |
Member Enumeration Documentation
◆ boundaryType
|
inherited |
Boundary type.
Either perfect magnetic conductor (PMC) or perfect electric conductor (PEC). Wether it's dirichlet or neumann boundary is dependent from the formulation.
Enumerator | |
---|---|
PMC | |
PEC | |
MAX_TYPE |
Definition at line 27 of file maxwell.hh.
◆ solverType
|
inherited |
Type of the solver.
Enumerator | |
---|---|
SUPERLU | |
SUPERLU2 | |
BICGSTAB | |
BICGSTAB2 | |
BICGSTABSUPERLU |
Definition at line 42 of file Maxwell2D_H.hh.
◆ subdivTypes
|
protectedinherited |
Definition at line 83 of file adaptiveModels.hh.
Constructor & Destructor Documentation
◆ Maxwell2D_H() [1/2]
hp2D::Maxwell2D_H::Maxwell2D_H | ( | concepts::EddyGeometry2D & | geom, |
enum concepts::MaxwellBoundary::boundaryType | bType = PMC , |
||
enum solverType | type = SUPERLU , |
||
bool | diagPrecond = true , |
||
bool | afterIter = false , |
||
const Real | eps = EPS0 , |
||
const Real | omega = OMEGA50 , |
||
const Real | mu = MU0 , |
||
const uint | geomRefAttrib = 100 |
||
) |
Constructor.
- Parameters
-
geom geometry, conductivity and source currents
bType type of boundary condition type solver type diagPrecond flag for preconditioning with diagonal matrix aterIter flag for after iterations eps dielectricity constant omega angular frequency in 1/s mu permeability constant in Ohm*s/m geomRefAttrib attrib for geometric refinement
◆ Maxwell2D_H() [2/2]
hp2D::Maxwell2D_H::Maxwell2D_H | ( | concepts::EddyGeometry2D & | geom, |
InputMaxwell2D_H & | input, | ||
const uint | geomRefAttrib = 100 |
||
) |
◆ ~Maxwell2D_H()
|
inlinevirtual |
Definition at line 144 of file Maxwell2D_H.hh.
Member Function Documentation
◆ bType() [1/2]
|
inlineinherited |
Definition at line 34 of file maxwell.hh.
◆ bType() [2/2]
|
inlineinherited |
Returns boundary type.
Definition at line 33 of file maxwell.hh.
◆ bTypeStr()
|
inlineinherited |
Returns name of boundary type as string.
Definition at line 36 of file maxwell.hh.
◆ constructor_()
|
privateinherited |
Private constructor.
◆ constructSpace_()
|
private |
Constructs the space.
◆ dissipation()
|
virtual |
Return dissipation power loss.
Implements concepts::MaxwellModel.
◆ eField()
concepts::ElementFormula<concepts::Cmplx2d>* hp2D::Maxwell2D_H::eField | ( | ) |
Returns a pointer to the e-Field.
Solves the problem, if not yet.
◆ hField()
concepts::ElementFormula<Cmplx>* hp2D::Maxwell2D_H::hField | ( | ) |
Returns a pointer to the h-Field (solution + h0).
Solves the problem, if not yet.
◆ identityMatrix_()
|
protectedinherited |
Calculate identity matrix and add's it to system matrix S
.
◆ info()
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from hp2D::Maxwell2D_H_Base.
◆ laplaceMatrix_()
|
protectedinherited |
Calculate stiffness matrix and add's it to system matrix S
.
◆ linearform_()
|
protectedinherited |
Calculate the load vector, assumes sigma = 0 inside the coil.
◆ magnEnergy()
|
virtual |
Return magnetic energy.
Implements concepts::MaxwellModel.
◆ matrices_()
|
private |
Building the matrices.
◆ mesh_()
|
inlineprotectedvirtualinherited |
◆ mshAbbr_()
|
inlineprotectedvirtualinherited |
Mesh abbreviation string.
Implements concepts::Model< Cmplx >.
Definition at line 67 of file Maxwell2D_H.hh.
◆ prebuild_()
|
inlineprotectedvirtual |
Space Prebuilder.
Implements hp2D::AdaptiveModel< Cmplx >.
Definition at line 162 of file Maxwell2D_H.hh.
◆ rebuildMesh() [1/2]
|
inherited |
◆ rebuildMesh() [2/2]
|
inherited |
Rebuilds only the mesh and sets the polynomial degrees.
- Parameters
-
l number of uniform refinements p number of polynomial enlargements g number of geometric refinements subdiv possibility to restrict subdivision strategy for geometric refinement
◆ rebuildMesh_()
|
protectedvirtualinherited |
Rebuilds only the mesh and sets the polynomial degrees.
Implements concepts::AdaptiveModel< Cmplx, 2 >.
◆ solution()
|
inlineinherited |
◆ solve_()
|
privatevirtual |
Method for solving, throws exception when it wasn't successfull.
Implements concepts::Model< Cmplx >.
◆ space() [1/2]
|
virtual |
Returns the space.
◆ space() [2/2]
|
pure virtualinherited |
Returns the space.
Friends And Related Function Documentation
◆ concepts::ModelControl< Maxwell2D_H >
|
friend |
Definition at line 112 of file Maxwell2D_H.hh.
Member Data Documentation
◆ A_
|
private |
Stiffness and System matrix.
Definition at line 167 of file Maxwell2D_H.hh.
◆ afterIter_
|
protectedinherited |
Nachiteration.
Definition at line 88 of file Maxwell2D_H.hh.
◆ bc_
|
protectedinherited |
Boundary conditions.
Definition at line 82 of file Maxwell2D_H.hh.
◆ bType_
|
privateinherited |
Type of boundary condition.
Definition at line 36 of file maxwell.hh.
◆ diagPrecond_
|
protectedinherited |
Using diagonal preconditioning.
Definition at line 86 of file Maxwell2D_H.hh.
◆ dissipation_
|
protectedinherited |
Dissipation power loss.
Definition at line 103 of file Maxwell2D_H.hh.
◆ eps_
|
protectedinherited |
Dielectricity constant.
Definition at line 97 of file Maxwell2D_H.hh.
◆ geom_
|
protectedinherited |
Mesh and material constants (sigma, j0)
Definition at line 80 of file Maxwell2D_H.hh.
◆ geomRefAttr_
|
protectedinherited |
Attribute of vertices or edges for geometric refinement.
Definition at line 91 of file adaptiveModels.hh.
◆ iOmegaEps_plus_Sigma_Inv_
|
protectedinherited |
◆ iterations_
|
protectedinherited |
Number of iterations for iterative solver.
Definition at line 109 of file Maxwell2D_H.hh.
◆ M_
|
private |
Mass matrix.
Definition at line 169 of file Maxwell2D_H.hh.
◆ magnEnergy_
|
protectedinherited |
Magnetic energy.
Definition at line 105 of file Maxwell2D_H.hh.
◆ matrixtime_
|
protectedinherited |
Definition at line 107 of file Maxwell2D_H.hh.
◆ mu_
|
protectedinherited |
Permeability constant.
Definition at line 101 of file Maxwell2D_H.hh.
◆ omega_
|
protectedinherited |
Angular frequency.
Definition at line 99 of file Maxwell2D_H.hh.
◆ problemName_
|
protectedinherited |
◆ residual_
|
protectedinherited |
Definition at line 89 of file Maxwell2D_H.hh.
◆ residualNorm_
|
protectedinherited |
Euclidian norm of the residual of solving the linear system.
Definition at line 91 of file Maxwell2D_H.hh.
◆ rhs_
|
protectedinherited |
Definition at line 93 of file Maxwell2D_H.hh.
◆ S_
|
private |
Definition at line 167 of file Maxwell2D_H.hh.
◆ sol_
|
protectedinherited |
◆ solvetime_
|
protectedinherited |
Time to solve the system, to build the matrices, to rebuild the space.
Definition at line 107 of file Maxwell2D_H.hh.
◆ spacetime_
|
protectedinherited |
Definition at line 107 of file Maxwell2D_H.hh.
◆ spc_
|
private |
Definition at line 165 of file Maxwell2D_H.hh.
◆ statusAfterIter_
|
protectedinherited |
Definition at line 88 of file Maxwell2D_H.hh.
◆ type_
|
protectedinherited |
Solver type.
Definition at line 82 of file Maxwell2D_H.hh.
The documentation for this class was generated from the following file:
- models/Maxwell2D_H.hh