12 #ifndef MFEM_JOULE_SOLVER
13 #define MFEM_JOULE_SOLVER
15 #include "../common/pfem_extras.hpp"
26 namespace electromagnetics
50 std::map<int, double> *materialMap;
60 if (materialMap != NULL) {
delete materialMap; }
187 std::map<int, double> sigmaAttMap,
188 std::map<int, double> TcapacityAttMap,
189 std::map<int, double> InvTcapAttMap,
190 std::map<int, double> InvTcondAttMap
220 void Debug(
const char *basefilename,
double time);
237 : E_gf(E_gf_), sigma(sigma_) {}
246 #endif // MFEM_USE_MPI
248 #endif // MFEM_JOULE_SOLVER
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
ParFiniteElementSpace & L2FESpace
void buildDiv(MeshDependentCoefficient &InvTcap)
void buildM2(MeshDependentCoefficient &alpha)
Class for grid function - Vector with associated FE space.
virtual ~MagneticDiffusionEOperator()
MeshDependentCoefficient * InvTcap
ParFiniteElementSpace & HDivFESpace
Base abstract class for first order time dependent operators.
Array< int > ess_bdr_vdofs
virtual ~MeshDependentCoefficient()
void e_exact(const Vector &x, double t, Vector &E)
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
void buildS1(double muInv)
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.
Abstract parallel finite element space.
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient in the element described by T at the point ip.
virtual ~ScaledGFCoefficient()
void Debug(const char *basefilename, double time)
void GetJouleHeating(ParGridFunction &E_gf, ParGridFunction &w_gf) const
void buildA1(double muInv, MeshDependentCoefficient &sigma, double dt)
ScaledGFCoefficient(GridFunction *gf, MeshDependentCoefficient &input_mdc)
double p_bc(const Vector &x, double t)
ParFiniteElementSpace & HCurlFESpace
void buildS2(MeshDependentCoefficient &alpha)
ParDiscreteLinearOperator * curl
Array< int > thermal_ess_bdr
void buildM3(MeshDependentCoefficient &Tcap)
void SetScaleFactor(const double &scale)
void buildM1(MeshDependentCoefficient &sigma)
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)
double muInv(const Vector &x)
ParMixedBilinearForm * weakDivC
MeshDependentCoefficient * Tcapacity
ParDiscreteLinearOperator * grad
double ElectricLosses(ParGridFunction &E_gf) const
JouleHeatingCoefficient(const MeshDependentCoefficient &sigma_, ParGridFunction &E_gf_)
MeshDependentCoefficient(const std::map< int, double > &inputMap, double scale=1.0)
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
void SetTime(const double t_)
Set the current time.
void buildA2(MeshDependentCoefficient &InvTcond, MeshDependentCoefficient &InvTcap, double dt)
MeshDependentCoefficient * sigma
Array< int > poisson_ess_bdr
Class for integration point with weight.
ParFiniteElementSpace & HGradFESpace
ParMixedBilinearForm * weakCurl
MeshDependentCoefficient * InvTcond
void SetMDC(const MeshDependentCoefficient &input_mdc)
Abstract class for hypre's solvers and preconditioners.
void edot_bc(const Vector &x, Vector &E)
double t_exact(const Vector &x)
Array< int > thermal_ess_bdr_vdofs
void b_exact(const Vector &x, double t, Vector &B)
Class for parallel grid function.
virtual ~JouleHeatingCoefficient()
virtual double Eval(ElementTransformation &T, const IntegrationPoint &ip)
Evaluate the coefficient at ip.
Wrapper for hypre's ParCSR matrix class.
void buildCurl(double muInv)
void buildA0(MeshDependentCoefficient &sigma)
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...
ParMixedBilinearForm * weakDiv
Array< int > poisson_ess_bdr_vdofs