12 #include "../../config/config.hpp"
21 using namespace miniapps;
23 namespace electromagnetics
26 VoltaSolver::VoltaSolver(
ParMesh & pmesh,
int order,
30 double (*phi_bc )(
const Vector&),
31 double (*rho_src)(
const Vector&),
68 MPI_Comm_size(pmesh_->
GetComm(), &num_procs_);
69 MPI_Comm_rank(pmesh_->
GetComm(), &myid_);
81 for (
int i=0; i<dbcs_->
Size(); i++)
83 ess_bdr_[(*dbcs_)[i]-1] = 1;
89 if ( phi_bc_ != NULL )
95 if ( rho_src_ != NULL )
101 if ( p_src_ != NULL )
144 if ( nbcs_->
Size() > 0 )
173 delete hCurlHDivEps_;
177 delete HCurlFESpace_;
179 map<string,socketstream*>::iterator mit;
180 for (mit=socks_.begin(); mit!=socks_.end(); mit++)
200 cout <<
"Number of H1 unknowns: " << size_h1 << endl;
201 cout <<
"Number of H(Curl) unknowns: " << size_nd << endl;
202 cout <<
"Number of H(Div) unknowns: " << size_rt << endl;
208 if (myid_ == 0) { cout <<
"Assembling ... " << flush; }
240 if (myid_ == 0) { cout <<
"done." << endl << flush; }
246 if (myid_ == 0) { cout <<
"Updating ..." << endl; }
251 H1FESpace_->
Update(
false);
252 HCurlFESpace_->
Update(
false);
253 HDivFESpace_->
Update(
false);
259 if ( rho_ ) { rho_->
Update(); }
260 if ( sigma_ ) { sigma_->
Update(); }
261 if ( p_ ) { p_->
Update(); }
268 if ( h1Mass_ ) { h1Mass_->
Update(); }
269 if ( h1SurfMass_ ) { h1SurfMass_->
Update(); }
270 if ( hCurlMass_ ) { hCurlMass_->
Update(); }
271 if ( hCurlHDiv_ ) { hCurlHDiv_->
Update(); }
280 if (myid_ == 0) { cout <<
"Running solver ... " << endl; }
285 if ( dbcs_->
Size() > 0 )
296 for (
int i=0; i<dbcs_->
Size(); i++)
300 dbc_bdr_attr[(*dbcs_)[i]-1] = 1;
318 MassH1->
Mult(*Rho,*RHS);
334 MassHCurl->
Mult(*P,*PD);
348 for (
int i=0; i<nbcs_->
Size(); i++)
352 nbc_bdr_attr[(*nbcs_)[i]-1] = 1;
359 MassS->
Mult(*Sigma,*RHS,1.0,1.0);
370 if ( dbcs_->
Size() > 0 )
397 pcg->
Mult(*RHS, *Phi);
413 Grad_->
Mult(*Phi,*E,-1.0);
419 if (myid_ == 0) { cout <<
"Computing D ..." << flush; }
422 hCurlHDivEps_->
Mult(*e_, ed);
425 hCurlHDiv_->
AddMult(*p_, ed, -1.0);
444 if (myid_ == 0) { cout <<
"done." << flush; }
450 if (myid_ == 0) { cout <<
"Solver done. " << endl; }
456 if (myid_ == 0) { cout <<
"Estimating Error ... " << flush; }
470 smooth_flux_fes, flux_fes, errors, norm_p);
472 if (myid_ == 0) { cout <<
"done." << endl; }
478 visit_dc_ = &visit_dc;
493 if (myid_ == 0) { cout <<
"Writing VisIt files ..." << flush; }
500 if (myid_ == 0) { cout <<
" done." << endl; }
507 if ( myid_ == 0 ) { cout <<
"Opening GLVis sockets." << endl; }
510 socks_[
"Phi"]->precision(8);
513 socks_[
"D"]->precision(8);
516 socks_[
"E"]->precision(8);
521 socks_[
"Rho"]->precision(8);
526 socks_[
"P"]->precision(8);
531 socks_[
"Sigma"]->precision(8);
538 if (myid_ == 0) { cout <<
"Sending data to GLVis ..." << flush; }
540 char vishost[] =
"localhost";
544 int Ww = 350, Wh = 350;
545 int offx = Ww+10, offy = Wh+45;
548 *phi_,
"Electric Potential (Phi)", Wx, Wy, Ww, Wh);
552 *d_,
"Electric Displacement (D)", Wx, Wy, Ww, Wh);
556 *e_,
"Electric Field (E)", Wx, Wy, Ww, Wh);
563 *rho_,
"Charge Density (Rho)", Wx, Wy, Ww, Wh);
569 *p_,
"Electric Polarization (P)", Wx, Wy, Ww, Wh);
575 *sigma_,
"Surface Charge Density (Sigma)", Wx, Wy, Ww, Wh);
578 if (myid_ == 0) { cout <<
" done." << endl; }
585 #endif // MFEM_USE_MPI
void WriteVisItFields(int it=0)
int Size() const
Logical size of the array.
void EliminateRowsCols(const Array< int > &rows_cols, const HypreParVector &X, HypreParVector &B)
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.
HYPRE_Int GetProblemSize()
virtual void Update(bool want_transform=true)
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 RegisterVisItFields(VisItDataCollection &visit_dc)
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)
void GetErrorEstimates(Vector &errors)
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.
Data collection with VisIt I/O routines.
void SetTime(double t)
Set physical time (for time-dependent simulations)
void SetMaxIter(int max_iter)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
int SpaceDimension() const
Wrapper for hypre's parallel vector class.
Array< int > bdr_attributes
Base class Coefficient that may optionally depend on time.
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
double L2ZZErrorEstimator(BilinearFormIntegrator &flux_integrator, const ParGridFunction &x, ParFiniteElementSpace &smooth_flux_fes, ParFiniteElementSpace &flux_fes, Vector &errors, int norm_p, double solver_tol, int solver_max_it)
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Abstract class for hypre's solvers and preconditioners.
class for C-function coefficient
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.
Integrator for (Q u, v) for VectorFiniteElements.
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Arbitrary order "L2-conforming" discontinuous finite elements.