12 #include "../../config/config.hpp"
23 using namespace miniapps;
25 namespace electromagnetics
28 TeslaSolver::TeslaSolver(
ParMesh & pmesh,
int order,
46 hDivHCurlMuInv_(NULL),
67 MPI_Comm_size(pmesh_->
GetComm(), &num_procs_);
68 MPI_Comm_rank(pmesh_->
GetComm(), &myid_);
83 for (
int i=0; i<kbcs.
Size(); i++)
85 non_k_bdr_[kbcs[i]-1] = 0;
104 if ( muInv_ == NULL )
114 if ( j_src_ != NULL )
121 if ( m_src_ != NULL )
155 if ( jCoef_ || kbcs.
Size() > 0 )
166 if ( kbcs.
Size() > 0 )
171 SurfCur_ =
new SurfaceCurrent(*H1FESpace_, *HCurlFESpace_, *Grad_,
207 delete curlMuInvCurl_;
209 delete hDivMassMuInv_;
210 delete hDivHCurlMuInv_;
213 delete HCurlFESpace_;
216 map<string,socketstream*>::iterator mit;
217 for (mit=socks_.begin(); mit!=socks_.end(); mit++)
237 cout <<
"Number of H1 unknowns: " << size_h1 << endl;
238 cout <<
"Number of H(Curl) unknowns: " << size_nd << endl;
239 cout <<
"Number of H(Div) unknowns: " << size_rt << endl << flush;
246 if (myid_ == 0) { cout <<
" Assembly ... " << flush; }
257 if ( j_ ) { j_->
Update(); }
258 if ( k_ ) { k_->
Update(); }
259 if ( m_ ) { m_->
Update(); }
273 if ( hDivMassMuInv_ )
280 if ( hDivHCurlMuInv_ )
282 hDivHCurlMuInv_->
Update();
289 if ( Grad_ ) { Grad_->
Update(); }
290 if ( DivFreeProj_ ) { DivFreeProj_->
Update(); }
291 if ( SurfCur_ ) { SurfCur_->
Update(); }
293 if (myid_ == 0) { cout <<
"done." << flush; }
299 if (myid_ == 0) { cout <<
"Running solver ... " << endl << flush; }
328 MassHCurl->
Mult(*J,*JD);
329 DivFreeProj_->
Mult(*JD, *RHS);
345 MassHDiv->
Mult(*M,*MD);
375 delete CurlMuInvCurl;
391 if (myid_ == 0) { cout <<
"Computing H ... " << flush; }
397 HDivHCurlMuInv->
Mult(*B,*BD);
401 HDivHCurlMuInv->
Mult(*M,*BD,-1.0*mu0_,1.0);
414 if (myid_ == 0) { cout <<
"done." << flush; }
418 delete HDivHCurlMuInv;
426 if (myid_ == 0) { cout <<
" Solver done. " << flush; }
432 if (myid_ == 0) { cout <<
"Estimating Error ... " << flush; }
445 smooth_flux_fes, flux_fes, errors, norm_p);
447 if (myid_ == 0) { cout <<
"done." << flush; }
453 visit_dc_ = &visit_dc;
469 if (myid_ == 0) { cout <<
"Writing VisIt files ..." << flush; }
476 if (myid_ == 0) { cout <<
" " << flush; }
483 if ( myid_ == 0 ) { cout <<
"Opening GLVis sockets." << endl << flush; }
486 socks_[
"A"]->precision(8);
489 socks_[
"B"]->precision(8);
492 socks_[
"H"]->precision(8);
497 socks_[
"J"]->precision(8);
502 socks_[
"K"]->precision(8);
505 socks_[
"Psi"]->precision(8);
510 socks_[
"M"]->precision(8);
512 if ( myid_ == 0 ) { cout <<
"GLVis sockets open." << endl << flush; }
518 if (myid_ == 0) { cout <<
"Sending data to GLVis ..." << flush; }
520 char vishost[] =
"localhost";
524 int Ww = 350, Wh = 350;
525 int offx = Ww+10, offy = Wh+45;
528 *a_,
"Vector Potential (A)", Wx, Wy, Ww, Wh);
532 *b_,
"Magnetic Flux Density (B)", Wx, Wy, Ww, Wh);
536 *h_,
"Magnetic Field (H)", Wx, Wy, Ww, Wh);
542 *j_,
"Current Density (J)", Wx, Wy, Ww, Wh);
550 *k_,
"Surface Current Density (K)", Wx, Wy, Ww, Wh);
555 "Surface Current Potential (Psi)", Wx, Wy, Ww, Wh);
561 *m_,
"Magnetization (M)", Wx, Wy, Ww, Wh);
564 if (myid_ == 0) { cout <<
" " << flush; }
572 : H1FESpace_(&H1FESpace),
573 HCurlFESpace_(&HCurlFESpace),
586 S0_ = s0_->ParallelAssemble();
599 for (
int i=0; i<vbcs_->
Size(); i++)
601 ess_bdr_[(*vbcs_)[i]-1] = 1;
606 for (
int i=0; i<kbcs_->
Size(); i++)
608 non_k_bdr_[(*kbcs_)[i]-1] = 0;
616 for (
int i=0; i<vbcs_->
Size(); i++)
620 vbc_bdr_attr[(*vbcs_)[i]-1] = 1;
628 s0_->ParallelEliminateEssentialBC(ess_bdr_,
648 if (myid_ == 0) { cout <<
"Computing K ... " << flush; }
651 pcg_->
Mult(*RHS_, *PSI_);
654 Grad_->
Mult(*PSI_,*K_);
657 Vector vZero(3); vZero = 0.0;
661 if (myid_ == 0) { cout <<
"done." << endl << flush; }
693 for (
int i=0; i<vbcs_->
Size(); i++)
697 vbc_bdr_attr[(*vbcs_)[i]-1] = 1;
714 #endif // MFEM_USE_MPI
int Size() const
Logical size of the array.
The Auxiliary-space Maxwell Solver in hypre.
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
virtual void RegisterField(const char *field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Subclass constant coefficient.
double muInv(const Vector &x)
Integrator for (curl u, curl v) for Nedelec elements.
HYPRE_Int GlobalTrueVSize()
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
virtual void Save()
Save the collection and a VisIt root file.
Abstract parallel finite element space.
void GetErrorEstimates(Vector &errors)
void ProjectCoefficient(Coefficient &coeff)
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
void SetPrintLevel(int print_lvl)
ParGridFunction * GetPsi()
The BoomerAMG solver in hypre.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
void VisualizeField(socketstream &sock, const char *vishost, int visport, ParGridFunction &gf, const char *title, int x, int y, int w, int h)
Jacobi preconditioner in hypre.
void SetPrintLevel(int print_level)
Data collection with VisIt I/O routines.
void L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, double solver_tol, int solver_max_it)
HYPRE_Int GetProblemSize()
void SetTime(double t)
Set physical time (for time-dependent simulations)
void RegisterVisItFields(VisItDataCollection &visit_dc)
void SetMaxIter(int max_iter)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
SurfaceCurrent(ParFiniteElementSpace &H1FESpace, ParFiniteElementSpace &HCurlFESpace, ParDiscreteGradOperator &Grad, Array< int > &kbcs, Array< int > &vbcs, Vector &vbcv)
int SpaceDimension() const
Wrapper for hypre's parallel vector class.
Array< int > bdr_attributes
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
void WriteVisItFields(int it=0)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
void ComputeSurfaceCurrent(ParGridFunction &k)
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Class for parallel grid function.
Wrapper for hypre's ParCSR matrix class.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Class for parallel meshes.
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Integrator for (Q u, v) for VectorFiniteElements.
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.