|
| 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, real_t mu, std::map< int, real_t > sigmaAttMap, std::map< int, real_t > TcapacityAttMap, std::map< int, real_t > InvTcapAttMap, std::map< int, real_t > 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.
|
|
virtual void | ImplicitSolve (const real_t dt, const Vector &x, Vector &k) |
| Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.
|
|
real_t | ElectricLosses (ParGridFunction &E_gf) const |
|
void | GetJouleHeating (ParGridFunction &E_gf, ParGridFunction &w_gf) const |
|
void | SetTime (const real_t t_) |
| Set the current time.
|
|
void | Debug (const char *basefilename, real_t time) |
|
virtual | ~MagneticDiffusionEOperator () |
|
| TimeDependentOperator (int n=0, real_t t_=0.0, Type type_=EXPLICIT) |
| Construct a "square" TimeDependentOperator y = f(x,t), where x and y have the same dimension n.
|
|
| TimeDependentOperator (int h, int w, real_t t_=0.0, Type type_=EXPLICIT) |
| Construct a TimeDependentOperator y = f(x,t), where x and y have dimensions w and h, respectively.
|
|
virtual real_t | GetTime () const |
| Read the currently set time.
|
|
bool | isExplicit () const |
| True if type is EXPLICIT.
|
|
bool | isImplicit () const |
| True if type is IMPLICIT or HOMOGENEOUS.
|
|
bool | isHomogeneous () const |
| True if type is HOMOGENEOUS.
|
|
EvalMode | GetEvalMode () const |
| Return the current evaluation mode. See SetEvalMode() for details.
|
|
virtual void | SetEvalMode (const EvalMode new_eval_mode) |
| Set the evaluation mode of the time-dependent operator.
|
|
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.
|
|
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.
|
|
virtual Operator & | GetImplicitGradient (const Vector &x, const Vector &k, real_t shift) const |
| Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time.
|
|
virtual Operator & | GetExplicitGradient (const Vector &x) const |
| Return an Operator representing dG/dx at the given point x and the currently set time.
|
|
virtual int | SUNImplicitSetup (const Vector &x, const Vector &fx, int jok, int *jcur, real_t gamma) |
| Setup the ODE linear system \( A(x,t) = (I - gamma J) \) or \( A = (M - gamma J) \), where \( J(x,t) = \frac{df}{dt(x,t)} \).
|
|
virtual int | SUNImplicitSolve (const Vector &b, Vector &x, real_t tol) |
| Solve the ODE linear system \( A x = b \) as setup by the method SUNImplicitSetup().
|
|
virtual int | SUNMassSetup () |
| Setup the mass matrix in the ODE system \( M y' = f(y,t) \) .
|
|
virtual int | SUNMassSolve (const Vector &b, Vector &x, real_t tol) |
| Solve the mass matrix linear system \( M x = b \) as setup by the method SUNMassSetup().
|
|
virtual int | SUNMassMult (const Vector &x, Vector &v) |
| Compute the mass matrix-vector product \( v = M x \) .
|
|
virtual | ~TimeDependentOperator () |
|
void | InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const |
| Initializes memory for true vectors of linear system.
|
|
| Operator (int s=0) |
| Construct a square Operator with given size s (default 0).
|
|
| Operator (int h, int w) |
| Construct an Operator with the given height (output size) and width (input size).
|
|
int | Height () const |
| Get the height (size of output) of the Operator. Synonym with NumRows().
|
|
int | NumRows () const |
| Get the number of rows (size of output) of the Operator. Synonym with Height().
|
|
int | Width () const |
| Get the width (size of input) of the Operator. Synonym with NumCols().
|
|
int | NumCols () const |
| Get the number of columns (size of input) of the Operator. Synonym with Width().
|
|
virtual MemoryClass | GetMemoryClass () const |
| Return the MemoryClass preferred by the Operator.
|
|
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.
|
|
virtual void | AddMult (const Vector &x, Vector &y, const real_t a=1.0) const |
| Operator application: y+=A(x) (default) or y+=a*A(x) .
|
|
virtual void | AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const |
| Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x) .
|
|
virtual void | ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const |
| Operator application on a matrix: Y=A(X) .
|
|
virtual void | ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const |
| Action of the transpose operator on a matrix: Y=A^t(X) .
|
|
virtual void | ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const |
| Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X) .
|
|
virtual void | ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const |
| Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X) .
|
|
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.
|
|
virtual void | AssembleDiagonal (Vector &diag) const |
| Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed.
|
|
virtual const Operator * | GetProlongation () const |
| Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity.
|
|
virtual const Operator * | GetRestriction () const |
| Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
|
|
virtual const Operator * | GetOutputProlongation () const |
| Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity.
|
|
virtual const Operator * | GetOutputRestrictionTranspose () const |
| Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators.
|
|
virtual const Operator * | GetOutputRestriction () const |
| Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
|
|
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.
|
|
void | FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B) |
| Form a column-constrained linear system using a matrix-free approach.
|
|
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() or Operator::FormRectangularLinearSystem().
|
|
void | FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A) |
| Return in A a parallel (on truedofs) version of this square operator.
|
|
void | FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A) |
| Return in A a parallel (on truedofs) version of this rectangular operator (including constraints).
|
|
void | FormDiscreteOperator (Operator *&A) |
| Return in A a parallel (on truedofs) version of this rectangular operator.
|
|
void | PrintMatlab (std::ostream &out, int n, int m=0) const |
| Prints operator with input size n and output size m in Matlab format.
|
|
virtual void | PrintMatlab (std::ostream &out) const |
| Prints operator in Matlab format.
|
|
virtual | ~Operator () |
| Virtual destructor.
|
|
Type | GetType () const |
| Return the type ID of the Operator class.
|
|
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 107 of file joule_solver.hpp.