MFEM
v3.3.2
Finite element discretization library
|
#include <joule_solver.hpp>
Public Member Functions | |
MagneticDiffusionEOperator (int len, ParFiniteElementSpace &L2FES, ParFiniteElementSpace &HCurlFES, ParFiniteElementSpace &HDivFES, ParFiniteElementSpace &HGradFES, Array< int > &ess_bdr, Array< int > &thermal_ess_bdr, Array< int > &poisson_ess_bdr, double mu, std::map< int, double > sigmaAttMap, std::map< int, double > TcapacityAttMap, std::map< int, double > InvTcapAttMap, std::map< int, double > InvTcondAttMap) | |
void | Init (Vector &vx) |
virtual void | Mult (const Vector &vx, Vector &dvx_dt) const |
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x, k, t) = G(x, t) and t is the current time. More... | |
virtual void | ImplicitSolve (const double dt, const Vector &x, Vector &k) |
Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t. More... | |
double | ElectricLosses (ParGridFunction &E_gf) const |
void | GetJouleHeating (ParGridFunction &E_gf, ParGridFunction &w_gf) const |
void | SetTime (const double _t) |
Set the current time. More... | |
void | Debug (const char *basefilename, double time) |
virtual | ~MagneticDiffusionEOperator () |
Public Member Functions inherited from mfem::TimeDependentOperator | |
TimeDependentOperator (int n=0, double t_=0.0, Type type_=EXPLICIT) | |
Construct a "square" TimeDependentOperator y = f(x,t), where x and y have the same dimension n. More... | |
TimeDependentOperator (int h, int w, double t_=0.0, Type type_=EXPLICIT) | |
Construct a TimeDependentOperator y = f(x,t), where x and y have dimensions w and h, respectively. More... | |
virtual double | GetTime () const |
Read the currently set time. More... | |
bool | isExplicit () const |
True if type is EXPLICIT. More... | |
bool | isImplicit () const |
True if type is IMPLICIT or HOMOGENEOUS. More... | |
bool | isHomogeneous () const |
True if type is HOMOGENEOUS. More... | |
virtual void | ExplicitMult (const Vector &x, Vector &y) const |
Perform the action of the explicit part of the operator, G: y = G(x, t) where t is the current time. More... | |
virtual void | ImplicitMult (const Vector &x, const Vector &k, Vector &y) const |
Perform the action of the implicit part of the operator, F: y = F(x, k, t) where t is the current time. More... | |
virtual Operator & | GetImplicitGradient (const Vector &x, const Vector &k, double shift) const |
Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time. More... | |
virtual Operator & | GetExplicitGradient (const Vector &x) const |
Return an Operator representing dG/dx at the given point x and the currently set time. More... | |
virtual | ~TimeDependentOperator () |
Public Member Functions inherited from mfem::Operator | |
Operator (int s=0) | |
Construct a square Operator with given size s (default 0). More... | |
Operator (int h, int w) | |
Construct an Operator with the given height (output size) and width (input size). More... | |
int | Height () const |
Get the height (size of output) of the Operator. Synonym with NumRows(). More... | |
int | NumRows () const |
Get the number of rows (size of output) of the Operator. Synonym with Height(). More... | |
int | Width () const |
Get the width (size of input) of the Operator. Synonym with NumCols(). More... | |
int | NumCols () const |
Get the number of columns (size of input) of the Operator. Synonym with Width(). More... | |
virtual void | MultTranspose (const Vector &x, Vector &y) const |
Action of the transpose operator: y=A^t(x) . The default behavior in class Operator is to generate an error. More... | |
virtual Operator & | GetGradient (const Vector &x) const |
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More... | |
virtual const Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More... | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More... | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0) |
Form a constrained linear system using a matrix-free approach. More... | |
virtual void | RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x) |
Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem(). More... | |
void | PrintMatlab (std::ostream &out, int n=0, int m=0) const |
Prints operator with input size n and output size m in Matlab format. More... | |
virtual | ~Operator () |
Virtual destructor. More... | |
Protected Member Functions | |
void | buildA0 (MeshDependentCoefficient &sigma) |
void | buildA1 (double muInv, MeshDependentCoefficient &sigma, double dt) |
void | buildA2 (MeshDependentCoefficient &InvTcond, MeshDependentCoefficient &InvTcap, double dt) |
void | buildM1 (MeshDependentCoefficient &sigma) |
void | buildM2 (MeshDependentCoefficient &alpha) |
void | buildM3 (MeshDependentCoefficient &Tcap) |
void | buildS1 (double muInv) |
void | buildS2 (MeshDependentCoefficient &alpha) |
void | buildGrad () |
void | buildCurl (double muInv) |
void | buildDiv (MeshDependentCoefficient &InvTcap) |
Additional Inherited Members | |
Public Types inherited from mfem::TimeDependentOperator | |
enum | Type { EXPLICIT, IMPLICIT, HOMOGENEOUS } |
Public Types inherited from mfem::Operator | |
enum | Type { MFEM_SPARSEMAT, Hypre_ParCSR, PETSC_MATAIJ, PETSC_MATIS, PETSC_MATSHELL, PETSC_MATNEST, PETSC_MATHYPRE, PETSC_MATGENERIC } |
Enumeration defining IDs for some classes derived from Operator. More... | |
After spatial discretization, the magnetic diffusion equation can be written as a system of ODEs:
S0(sigma) P = 0 dE/dt = - (M1(sigma) + dt S1(1/mu))^{-1}*(S1(1/mu)*E + sigma Grad P) dB/dt = - Curl(E) dF/dt = (M2(c/k) + dt S2(1))^{-1} (-S2(1) F + Div J) dcT/dt = -Div F + J
where
P is the 0-form electrostatic potential, E is the 1-form electric field, B is the 2-form magnetic flux, F is the 2-form thermal flux, T is the 3-form temperature. M is the mass matrix, S is the stiffness matrix, Curl is the curl matrix, Div is the divergence matrix, sigma is the electrical conductivity, k is the thermal conductivity, c is the heat capacity, J is a function of the Joule heating sigma (E dot E).
Class MagneticDiffusionEOperator represents the right-hand side of the above system of ODEs.
Definition at line 106 of file joule_solver.hpp.
mfem::electromagnetics::MagneticDiffusionEOperator::MagneticDiffusionEOperator | ( | int | len, |
ParFiniteElementSpace & | L2FES, | ||
ParFiniteElementSpace & | HCurlFES, | ||
ParFiniteElementSpace & | HDivFES, | ||
ParFiniteElementSpace & | HGradFES, | ||
Array< int > & | ess_bdr, | ||
Array< int > & | thermal_ess_bdr, | ||
Array< int > & | poisson_ess_bdr, | ||
double | mu, | ||
std::map< int, double > | sigmaAttMap, | ||
std::map< int, double > | TcapacityAttMap, | ||
std::map< int, double > | InvTcapAttMap, | ||
std::map< int, double > | InvTcondAttMap | ||
) |
Definition at line 26 of file joule_solver.cpp.
|
virtual |
Definition at line 848 of file joule_solver.cpp.
|
protected |
Definition at line 664 of file joule_solver.cpp.
|
protected |
Definition at line 681 of file joule_solver.cpp.
|
protected |
Definition at line 703 of file joule_solver.cpp.
|
protected |
Definition at line 776 of file joule_solver.cpp.
|
protected |
Definition at line 793 of file joule_solver.cpp.
|
protected |
Definition at line 809 of file joule_solver.cpp.
|
protected |
Definition at line 721 of file joule_solver.cpp.
|
protected |
Definition at line 732 of file joule_solver.cpp.
|
protected |
Definition at line 744 of file joule_solver.cpp.
|
protected |
Definition at line 756 of file joule_solver.cpp.
|
protected |
Definition at line 766 of file joule_solver.cpp.
void mfem::electromagnetics::MagneticDiffusionEOperator::Debug | ( | const char * | basefilename, |
double | time | ||
) |
Definition at line 894 of file joule_solver.cpp.
double mfem::electromagnetics::MagneticDiffusionEOperator::ElectricLosses | ( | ParGridFunction & | E_gf | ) | const |
Definition at line 820 of file joule_solver.cpp.
void mfem::electromagnetics::MagneticDiffusionEOperator::GetJouleHeating | ( | ParGridFunction & | E_gf, |
ParGridFunction & | w_gf | ||
) | const |
Definition at line 833 of file joule_solver.cpp.
|
virtual |
Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.
For general F and G, the equation for k becomes: F(x + dt k, k, t) = G(x + dt k, t).
The input vector x corresponds to time index (or cycle) n, while the currently set time, t, and the result vector k correspond to time index n+1. The time step dt corresponds to the time interval between cycles n and n+1.
This method allows for the abstract implementation of some time integration methods, including diagonal implicit Runge-Kutta (DIRK) methods and the backward Euler method in particular.
If not re-implemented, this method simply generates an error.
Reimplemented from mfem::TimeDependentOperator.
Definition at line 412 of file joule_solver.cpp.
void mfem::electromagnetics::MagneticDiffusionEOperator::Init | ( | Vector & | vx | ) |
Definition at line 101 of file joule_solver.cpp.
|
virtual |
Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x, k, t) = G(x, t) and t is the current time.
Reimplemented from mfem::TimeDependentOperator.
Definition at line 162 of file joule_solver.cpp.
|
virtual |
Set the current time.
Reimplemented from mfem::TimeDependentOperator.
Definition at line 845 of file joule_solver.cpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 127 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 127 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 127 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 139 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 137 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 141 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.
|
protected |
Definition at line 122 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 145 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 147 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 143 of file joule_solver.hpp.
|
protected |
Definition at line 158 of file joule_solver.hpp.
|
protected |
Definition at line 158 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 150 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 151 of file joule_solver.hpp.
|
protected |
Definition at line 122 of file joule_solver.hpp.
|
protected |
Definition at line 115 of file joule_solver.hpp.
|
protected |
Definition at line 116 of file joule_solver.hpp.
|
protected |
Definition at line 117 of file joule_solver.hpp.
|
protected |
Definition at line 157 of file joule_solver.hpp.
|
protected |
Definition at line 157 of file joule_solver.hpp.
|
protected |
Definition at line 114 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 127 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 127 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 127 of file joule_solver.hpp.
|
protected |
Definition at line 158 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 138 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 142 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 140 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 146 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 148 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 144 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 154 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 155 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 121 of file joule_solver.hpp.
|
protected |
Definition at line 157 of file joule_solver.hpp.
|
protected |
Definition at line 157 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 152 of file joule_solver.hpp.
|
mutableprotected |
Definition at line 153 of file joule_solver.hpp.
|
protected |
Definition at line 131 of file joule_solver.hpp.
|
protected |
Definition at line 131 of file joule_solver.hpp.
|
protected |
Definition at line 131 of file joule_solver.hpp.
|
protected |
Definition at line 123 of file joule_solver.hpp.
|
protected |
Definition at line 123 of file joule_solver.hpp.
|
protected |
Definition at line 123 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.
|
protected |
Definition at line 128 of file joule_solver.hpp.