MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tesla_solver.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "../../config/config.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 #include "tesla_solver.hpp"
17 
18 using namespace std;
19 
20 namespace mfem
21 {
22 
23 using namespace miniapps;
24 
25 namespace electromagnetics
26 {
27 
28 TeslaSolver::TeslaSolver(ParMesh & pmesh, int order,
29  Array<int> & kbcs,
30  Array<int> & vbcs, Vector & vbcv,
31  double (*muInv)(const Vector&),
32  void (*a_bc )(const Vector&, Vector&),
33  void (*j_src)(const Vector&, Vector&),
34  void (*m_src)(const Vector&, Vector&))
35  : myid_(0),
36  num_procs_(1),
37  order_(order),
38  pmesh_(&pmesh),
39  visit_dc_(NULL),
40  H1FESpace_(NULL),
41  HCurlFESpace_(NULL),
42  HDivFESpace_(NULL),
43  curlMuInvCurl_(NULL),
44  hCurlMass_(NULL),
45  hDivMassMuInv_(NULL),
46  hDivHCurlMuInv_(NULL),
47  Grad_(NULL),
48  Curl_(NULL),
49  a_(NULL),
50  b_(NULL),
51  h_(NULL),
52  j_(NULL),
53  k_(NULL),
54  m_(NULL),
55  DivFreeProj_(NULL),
56  SurfCur_(NULL),
57  muInvCoef_(NULL),
58  aBCCoef_(NULL),
59  jCoef_(NULL),
60  mCoef_(NULL),
61  muInv_(muInv),
62  a_bc_(a_bc),
63  j_src_(j_src),
64  m_src_(m_src)
65 {
66  // Initialize MPI variables
67  MPI_Comm_size(pmesh_->GetComm(), &num_procs_);
68  MPI_Comm_rank(pmesh_->GetComm(), &myid_);
69 
70  // Define compatible parallel finite element spaces on the parallel
71  // mesh. Here we use arbitrary order H1, Nedelec, and Raviart-Thomas finite
72  // elements.
73  H1FESpace_ = new H1_ParFESpace(pmesh_,order,pmesh_->Dimension());
74  HCurlFESpace_ = new ND_ParFESpace(pmesh_,order,pmesh_->Dimension());
75  HDivFESpace_ = new RT_ParFESpace(pmesh_,order,pmesh_->Dimension());
76 
77  // Select surface attributes for Dirichlet BCs
78  ess_bdr_.SetSize(pmesh.bdr_attributes.Max());
79  non_k_bdr_.SetSize(pmesh.bdr_attributes.Max());
80  ess_bdr_ = 1; // All outer surfaces
81  non_k_bdr_ = 1; // Surfaces without applied surface currents
82 
83  for (int i=0; i<kbcs.Size(); i++)
84  {
85  non_k_bdr_[kbcs[i]-1] = 0;
86  }
87 
88  // Setup various coefficients
89 
90  // Vector Potential on the outer surface
91  if ( a_bc_ == NULL )
92  {
93  Vector Zero(3);
94  Zero = 0.0;
95  aBCCoef_ = new VectorConstantCoefficient(Zero);
96  }
97  else
98  {
99  aBCCoef_ = new VectorFunctionCoefficient(pmesh_->SpaceDimension(),
100  *a_bc_);
101  }
102 
103  // Inverse of the magnetic permeability
104  if ( muInv_ == NULL )
105  {
106  muInvCoef_ = new ConstantCoefficient(1.0/mu0_);
107  }
108  else
109  {
110  muInvCoef_ = new FunctionCoefficient(muInv_);
111  }
112 
113  // Volume Current Density
114  if ( j_src_ != NULL )
115  {
116  jCoef_ = new VectorFunctionCoefficient(pmesh_->SpaceDimension(),
117  j_src_);
118  }
119 
120  // Magnetization
121  if ( m_src_ != NULL )
122  {
123  mCoef_ = new VectorFunctionCoefficient(pmesh_->SpaceDimension(),
124  m_src_);
125  }
126 
127  // Bilinear Forms
128  curlMuInvCurl_ = new ParBilinearForm(HCurlFESpace_);
129  curlMuInvCurl_->AddDomainIntegrator(new CurlCurlIntegrator(*muInvCoef_));
130 
131  hCurlMass_ = new ParBilinearForm(HCurlFESpace_);
133 
134  hDivHCurlMuInv_ = new ParMixedBilinearForm(HDivFESpace_, HCurlFESpace_);
135  hDivHCurlMuInv_->AddDomainIntegrator(new VectorFEMassIntegrator(*muInvCoef_));
136 
137  // Assemble Matrices
138  curlMuInvCurl_->Assemble();
139  curlMuInvCurl_->Finalize();
140 
141  hCurlMass_->Assemble();
142  hCurlMass_->Finalize();
143 
144  hDivHCurlMuInv_->Assemble();
145  hDivHCurlMuInv_->Finalize();
146 
147  // Discrete Curl operator
148  Curl_ = new ParDiscreteCurlOperator(HCurlFESpace_, HDivFESpace_);
149 
150  // Build grid functions
151  a_ = new ParGridFunction(HCurlFESpace_);
152  b_ = new ParGridFunction(HDivFESpace_);
153  h_ = new ParGridFunction(HCurlFESpace_);
154 
155  if ( jCoef_ || kbcs.Size() > 0 )
156  {
157  Grad_ = new ParDiscreteGradOperator(H1FESpace_, HCurlFESpace_);
158  }
159  if ( jCoef_ )
160  {
161  j_ = new ParGridFunction(HCurlFESpace_);
162  DivFreeProj_ = new DivergenceFreeProjector(*H1FESpace_, *HCurlFESpace_,
163  *Grad_);
164  }
165 
166  if ( kbcs.Size() > 0 )
167  {
168  k_ = new ParGridFunction(HCurlFESpace_);
169 
170  // Object to solve the subproblem of computing surface currents
171  SurfCur_ = new SurfaceCurrent(*H1FESpace_, *HCurlFESpace_, *Grad_,
172  kbcs, vbcs, vbcv);
173  }
174 
175  if ( mCoef_ )
176  {
177  m_ = new ParGridFunction(HDivFESpace_);
178 
179  hDivMassMuInv_ = new ParBilinearForm(HDivFESpace_);
180  hDivMassMuInv_->AddDomainIntegrator(
181  new VectorFEMassIntegrator(*muInvCoef_));
182  hDivMassMuInv_->Assemble();
183  hDivMassMuInv_->Finalize();
184  }
185 }
186 
188 {
189  delete muInvCoef_;
190  delete jCoef_;
191  delete mCoef_;
192  delete aBCCoef_;
193 
194  delete DivFreeProj_;
195  delete SurfCur_;
196 
197  delete a_;
198  delete b_;
199  delete h_;
200  delete j_;
201  delete k_;
202  delete m_;
203 
204  delete Grad_;
205  delete Curl_;
206 
207  delete curlMuInvCurl_;
208  delete hCurlMass_;
209  delete hDivMassMuInv_;
210  delete hDivHCurlMuInv_;
211 
212  delete H1FESpace_;
213  delete HCurlFESpace_;
214  delete HDivFESpace_;
215 
216  map<string,socketstream*>::iterator mit;
217  for (mit=socks_.begin(); mit!=socks_.end(); mit++)
218  {
219  delete mit->second;
220  }
221 }
222 
223 HYPRE_Int
225 {
226  return HCurlFESpace_->GlobalTrueVSize();
227 }
228 
229 void
231 {
232  HYPRE_Int size_h1 = H1FESpace_->GlobalTrueVSize();
233  HYPRE_Int size_nd = HCurlFESpace_->GlobalTrueVSize();
234  HYPRE_Int size_rt = HDivFESpace_->GlobalTrueVSize();
235  if (myid_ == 0)
236  {
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;
240  }
241 }
242 
243 void
245 {
246  if (myid_ == 0) { cout << " Assembly ... " << flush; }
247 
248  // Inform the spaces that the mesh has changed
249  H1FESpace_->Update();
250  HCurlFESpace_->Update();
251  HDivFESpace_->Update();
252 
253  // Inform the grid functions that the space has changed.
254  a_->Update();
255  h_->Update();
256  b_->Update();
257  if ( j_ ) { j_->Update(); }
258  if ( k_ ) { k_->Update(); }
259  if ( m_ ) { m_->Update(); }
260 
261  // Inform the bilinear forms that the space has changed.
262  curlMuInvCurl_->Update();
263  curlMuInvCurl_->Assemble();
264  curlMuInvCurl_->Finalize();
265 
266  if ( hCurlMass_ )
267  {
268  hCurlMass_->Update();
269  hCurlMass_->Assemble();
270  hCurlMass_->Finalize();
271  }
272 
273  if ( hDivMassMuInv_ )
274  {
275  hDivMassMuInv_->Update();
276  hDivMassMuInv_->Assemble();
277  hDivMassMuInv_->Finalize();
278  }
279 
280  if ( hDivHCurlMuInv_ )
281  {
282  hDivHCurlMuInv_->Update();
283  hDivHCurlMuInv_->Assemble();
284  hDivHCurlMuInv_->Finalize();
285  }
286 
287  // Inform the other objects that the space has changed.
288  Curl_->Update();
289  if ( Grad_ ) { Grad_->Update(); }
290  if ( DivFreeProj_ ) { DivFreeProj_->Update(); }
291  if ( SurfCur_ ) { SurfCur_->Update(); }
292 
293  if (myid_ == 0) { cout << "done." << flush; }
294 }
295 
296 void
298 {
299  if (myid_ == 0) { cout << "Running solver ... " << endl << flush; }
300 
301  // Initialize the magnetic vector potential with its boundary conditions
302  *a_ = 0.0;
303 
304  // Apply surface currents if available
305  if ( k_ )
306  {
307  SurfCur_->ComputeSurfaceCurrent(*k_);
308  *a_ = *k_;
309  }
310 
311  // Apply uniform B boundary condition on remaining surfaces
312  a_->ProjectBdrCoefficientTangent(*aBCCoef_, non_k_bdr_);
313 
314  // Initialize the RHS vector
315  HypreParVector *RHS = new HypreParVector(HCurlFESpace_);
316  *RHS = 0.0;
317 
318  HypreParMatrix *MassHCurl = hCurlMass_->ParallelAssemble();
319 
320  // Initialize the volumetric current density
321  if ( j_ )
322  {
323  j_->ProjectCoefficient(*jCoef_);
324 
325  HypreParVector *J = j_->ParallelProject();
326  HypreParVector *JD = new HypreParVector(HCurlFESpace_);
327 
328  MassHCurl->Mult(*J,*JD);
329  DivFreeProj_->Mult(*JD, *RHS);
330 
331  delete J;
332  delete JD;
333  }
334 
335  // Initialize the Magnetization
336  HypreParVector *M = NULL;
337  if ( m_ )
338  {
339  m_->ProjectCoefficient(*mCoef_);
340  M = m_->ParallelProject();
341 
342  HypreParMatrix *MassHDiv = hDivMassMuInv_->ParallelAssemble();
343  HypreParVector *MD = new HypreParVector(HDivFESpace_);
344 
345  MassHDiv->Mult(*M,*MD);
346  Curl_->MultTranspose(*MD,*RHS,mu0_,1.0);
347 
348  delete MassHDiv;
349  delete MD;
350  }
351 
352  // Apply Dirichlet BCs to matrix and right hand side
353  HypreParMatrix *CurlMuInvCurl = curlMuInvCurl_->ParallelAssemble();
354  HypreParVector *A = a_->ParallelProject();
355 
356  // Apply the boundary conditions to the assembled matrix and vectors
357  curlMuInvCurl_->ParallelEliminateEssentialBC(ess_bdr_,
358  *CurlMuInvCurl,
359  *A, *RHS);
360 
361  // Define and apply a parallel PCG solver for AX=B with the AMS
362  // preconditioner from hypre.
363  HypreAMS *ams = new HypreAMS(*CurlMuInvCurl, HCurlFESpace_);
364  ams->SetSingularProblem();
365 
366  HyprePCG *pcg = new HyprePCG(*CurlMuInvCurl);
367  pcg->SetTol(1e-12);
368  pcg->SetMaxIter(500);
369  pcg->SetPrintLevel(2);
370  pcg->SetPreconditioner(*ams);
371  pcg->Mult(*RHS, *A);
372 
373  delete ams;
374  delete pcg;
375  delete CurlMuInvCurl;
376  delete RHS;
377 
378  // Extract the parallel grid function corresponding to the finite
379  // element approximation Phi. This is the local solution on each
380  // processor.
381  *a_ = *A;
382 
383  // Compute the negative Gradient of the solution vector. This is
384  // the magnetic field corresponding to the scalar potential
385  // represented by phi.
386  HypreParVector *B = new HypreParVector(HDivFESpace_);
387  Curl_->Mult(*A,*B);
388  *b_ = *B;
389 
390  // Compute magnetic field (H) from B and M
391  if (myid_ == 0) { cout << "Computing H ... " << flush; }
392 
393  HypreParMatrix *HDivHCurlMuInv = hDivHCurlMuInv_->ParallelAssemble();
394  HypreParVector *BD = new HypreParVector(HCurlFESpace_);
395  HypreParVector *H = new HypreParVector(HCurlFESpace_);
396 
397  HDivHCurlMuInv->Mult(*B,*BD);
398 
399  if ( M )
400  {
401  HDivHCurlMuInv->Mult(*M,*BD,-1.0*mu0_,1.0);
402  }
403 
404  HyprePCG * pcgM = new HyprePCG(*MassHCurl);
405  pcgM->SetTol(1e-12);
406  pcgM->SetMaxIter(500);
407  pcgM->SetPrintLevel(0);
408  HypreDiagScale *diagM = new HypreDiagScale;
409  pcgM->SetPreconditioner(*diagM);
410  pcgM->Mult(*BD,*H);
411 
412  *h_ = *H;
413 
414  if (myid_ == 0) { cout << "done." << flush; }
415 
416  delete diagM;
417  delete pcgM;
418  delete HDivHCurlMuInv;
419  delete MassHCurl;
420  delete A;
421  delete B;
422  delete BD;
423  delete H;
424  delete M;
425 
426  if (myid_ == 0) { cout << " Solver done. " << flush; }
427 }
428 
429 void
431 {
432  if (myid_ == 0) { cout << "Estimating Error ... " << flush; }
433 
434  // Space for the discontinuous (original) flux
435  CurlCurlIntegrator flux_integrator(*muInvCoef_);
436  RT_FECollection flux_fec(order_-1, pmesh_->SpaceDimension());
437  ParFiniteElementSpace flux_fes(pmesh_, &flux_fec);
438 
439  // Space for the smoothed (conforming) flux
440  double norm_p = 1;
441  ND_FECollection smooth_flux_fec(order_, pmesh_->Dimension());
442  ParFiniteElementSpace smooth_flux_fes(pmesh_, &smooth_flux_fec);
443 
444  L2ZZErrorEstimator(flux_integrator, *a_,
445  smooth_flux_fes, flux_fes, errors, norm_p);
446 
447  if (myid_ == 0) { cout << "done." << flush; }
448 }
449 
450 void
452 {
453  visit_dc_ = &visit_dc;
454 
455  visit_dc.RegisterField("A", a_);
456  visit_dc.RegisterField("B", b_);
457  visit_dc.RegisterField("H", h_);
458  if ( j_ ) { visit_dc.RegisterField("J", j_); }
459  if ( k_ ) { visit_dc.RegisterField("K", k_); }
460  if ( m_ ) { visit_dc.RegisterField("M", m_); }
461  if ( SurfCur_ ) { visit_dc.RegisterField("Psi", SurfCur_->GetPsi()); }
462 }
463 
464 void
466 {
467  if ( visit_dc_ )
468  {
469  if (myid_ == 0) { cout << "Writing VisIt files ..." << flush; }
470 
471  HYPRE_Int prob_size = this->GetProblemSize();
472  visit_dc_->SetCycle(it);
473  visit_dc_->SetTime(prob_size);
474  visit_dc_->Save();
475 
476  if (myid_ == 0) { cout << " " << flush; }
477  }
478 }
479 
480 void
482 {
483  if ( myid_ == 0 ) { cout << "Opening GLVis sockets." << endl << flush; }
484 
485  socks_["A"] = new socketstream;
486  socks_["A"]->precision(8);
487 
488  socks_["B"] = new socketstream;
489  socks_["B"]->precision(8);
490 
491  socks_["H"] = new socketstream;
492  socks_["H"]->precision(8);
493 
494  if ( j_ )
495  {
496  socks_["J"] = new socketstream;
497  socks_["J"]->precision(8);
498  }
499  if ( k_ )
500  {
501  socks_["K"] = new socketstream;
502  socks_["K"]->precision(8);
503 
504  socks_["Psi"] = new socketstream;
505  socks_["Psi"]->precision(8);
506  }
507  if ( m_ )
508  {
509  socks_["M"] = new socketstream;
510  socks_["M"]->precision(8);
511  }
512  if ( myid_ == 0 ) { cout << "GLVis sockets open." << endl << flush; }
513 }
514 
515 void
517 {
518  if (myid_ == 0) { cout << "Sending data to GLVis ..." << flush; }
519 
520  char vishost[] = "localhost";
521  int visport = 19916;
522 
523  int Wx = 0, Wy = 0; // window position
524  int Ww = 350, Wh = 350; // window size
525  int offx = Ww+10, offy = Wh+45; // window offsets
526 
527  VisualizeField(*socks_["A"], vishost, visport,
528  *a_, "Vector Potential (A)", Wx, Wy, Ww, Wh);
529  Wx += offx;
530 
531  VisualizeField(*socks_["B"], vishost, visport,
532  *b_, "Magnetic Flux Density (B)", Wx, Wy, Ww, Wh);
533  Wx += offx;
534 
535  VisualizeField(*socks_["H"], vishost, visport,
536  *h_, "Magnetic Field (H)", Wx, Wy, Ww, Wh);
537  Wx += offx;
538 
539  if ( j_ )
540  {
541  VisualizeField(*socks_["J"], vishost, visport,
542  *j_, "Current Density (J)", Wx, Wy, Ww, Wh);
543  }
544 
545  Wx = 0; Wy += offy; // next line
546 
547  if ( k_ )
548  {
549  VisualizeField(*socks_["K"], vishost, visport,
550  *k_, "Surface Current Density (K)", Wx, Wy, Ww, Wh);
551  Wx += offx;
552 
553  VisualizeField(*socks_["Psi"], vishost, visport,
554  *SurfCur_->GetPsi(),
555  "Surface Current Potential (Psi)", Wx, Wy, Ww, Wh);
556  Wx += offx;
557  }
558  if ( m_ )
559  {
560  VisualizeField(*socks_["M"], vishost, visport,
561  *m_, "Magnetization (M)", Wx, Wy, Ww, Wh);
562  Wx += offx;
563  }
564  if (myid_ == 0) { cout << " " << flush; }
565 }
566 
568  ParFiniteElementSpace & HCurlFESpace,
570  Array<int> & kbcs,
571  Array<int> & vbcs, Vector & vbcv)
572  : H1FESpace_(&H1FESpace),
573  HCurlFESpace_(&HCurlFESpace),
574  Grad_(&Grad),
575  kbcs_(&kbcs),
576  vbcs_(&vbcs),
577  vbcv_(&vbcv)
578 {
579  // Initialize MPI variables
580  MPI_Comm_rank(H1FESpace_->GetParMesh()->GetComm(), &myid_);
581 
582  s0_ = new ParBilinearForm(H1FESpace_);
583  s0_->AddBoundaryIntegrator(new DiffusionIntegrator);
584  s0_->Assemble();
585  s0_->Finalize();
586  S0_ = s0_->ParallelAssemble();
587 
588  amg_ = new HypreBoomerAMG(*S0_);
589  amg_->SetPrintLevel(0);
590 
591  pcg_ = new HyprePCG(*S0_);
592  pcg_->SetTol(1e-12);
593  pcg_->SetMaxIter(200);
594  pcg_->SetPrintLevel(0);
595  pcg_->SetPreconditioner(*amg_);
596 
597  ess_bdr_.SetSize(H1FESpace_->GetParMesh()->bdr_attributes.Max());
598  ess_bdr_ = 0;
599  for (int i=0; i<vbcs_->Size(); i++)
600  {
601  ess_bdr_[(*vbcs_)[i]-1] = 1;
602  }
603 
604  non_k_bdr_.SetSize(H1FESpace_->GetParMesh()->bdr_attributes.Max());
605  non_k_bdr_ = 1;
606  for (int i=0; i<kbcs_->Size(); i++)
607  {
608  non_k_bdr_[(*kbcs_)[i]-1] = 0;
609  }
610 
611  psi_ = new ParGridFunction(H1FESpace_);
612  *psi_ = 0.0;
613 
614  // Apply piecewise constant voltage boundary condition
615  Array<int> vbc_bdr_attr(H1FESpace_->GetParMesh()->bdr_attributes.Max());
616  for (int i=0; i<vbcs_->Size(); i++)
617  {
618  ConstantCoefficient voltage((*vbcv_)[i]);
619  vbc_bdr_attr = 0;
620  vbc_bdr_attr[(*vbcs_)[i]-1] = 1;
621  psi_->ProjectBdrCoefficient(voltage, vbc_bdr_attr);
622  }
623 
624  PSI_ = psi_->ParallelProject();
625  RHS_ = new HypreParVector(H1FESpace_);
626  K_ = new HypreParVector(HCurlFESpace_);
627 
628  s0_->ParallelEliminateEssentialBC(ess_bdr_,
629  *S0_,
630  *PSI_, *RHS_);
631 }
632 
634 {
635  delete pcg_;
636  delete amg_;
637  delete S0_;
638  delete PSI_;
639  delete RHS_;
640  delete K_;
641  delete s0_;
642  delete psi_;
643 }
644 
645 void
647 {
648  if (myid_ == 0) { cout << "Computing K ... " << flush; }
649 
650  k = 0.0;
651  pcg_->Mult(*RHS_, *PSI_);
652  *psi_ = *PSI_;
653 
654  Grad_->Mult(*PSI_,*K_);
655  k = *K_;
656 
657  Vector vZero(3); vZero = 0.0;
658  VectorConstantCoefficient Zero(vZero);
659  k.ProjectBdrCoefficientTangent(Zero,non_k_bdr_);
660 
661  if (myid_ == 0) { cout << "done." << endl << flush; }
662 }
663 
664 void
666 {
667  delete pcg_;
668  delete amg_;
669  delete S0_;
670  delete PSI_;
671  delete RHS_;
672  delete K_;
673 
674  psi_->Update();
675  *psi_ = 0.0;
676 
677  s0_->Update();
678  s0_->Assemble();
679  s0_->Finalize();
680  S0_ = s0_->ParallelAssemble();
681 
682  amg_ = new HypreBoomerAMG(*S0_);
683  amg_->SetPrintLevel(0);
684 
685  pcg_ = new HyprePCG(*S0_);
686  pcg_->SetTol(1e-12);
687  pcg_->SetMaxIter(200);
688  pcg_->SetPrintLevel(0);
689  pcg_->SetPreconditioner(*amg_);
690 
691  // Apply piecewise constant voltage boundary condition
692  Array<int> vbc_bdr_attr(H1FESpace_->GetParMesh()->bdr_attributes.Max());
693  for (int i=0; i<vbcs_->Size(); i++)
694  {
695  ConstantCoefficient voltage((*vbcv_)[i]);
696  vbc_bdr_attr = 0;
697  vbc_bdr_attr[(*vbcs_)[i]-1] = 1;
698  psi_->ProjectBdrCoefficient(voltage, vbc_bdr_attr);
699  }
700 
701  PSI_ = psi_->ParallelProject();
702  RHS_ = new HypreParVector(H1FESpace_);
703  K_ = new HypreParVector(HCurlFESpace_);
704 
705  s0_->ParallelEliminateEssentialBC(ess_bdr_,
706  *S0_,
707  *PSI_, *RHS_);
708 }
709 
710 } // namespace electromagnetics
711 
712 } // namespace mfem
713 
714 #endif // MFEM_USE_MPI
void SetTol(double tol)
Definition: hypre.cpp:1917
int Size() const
Logical size of the array.
Definition: array.hpp:109
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:761
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.
Definition: coefficient.hpp:57
double muInv(const Vector &x)
Definition: tesla.cpp:70
Integrator for (curl u, curl v) for Nedelec elements.
Definition: bilininteg.hpp:401
HYPRE_Int MultTranspose(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A^t * x + beta * y.
Definition: pfem_extras.cpp:98
virtual void Save()
Save the collection and a VisIt root file.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void GetErrorEstimates(Vector &errors)
void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:241
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:903
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:1932
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
virtual void Update(FiniteElementSpace *nfes=NULL)
HypreParMatrix * ParallelAssemble()
Returns the matrix assembled on the true dofs, i.e. P^t A P.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
The BoomerAMG solver in hypre.
Definition: hypre.hpp:698
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0)
Computes y = alpha * A * x + beta * y.
Definition: pfem_extras.cpp:82
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.
Definition: hypre.hpp:659
void SetPrintLevel(int print_level)
Definition: hypre.hpp:739
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)
Definition: pgridfunc.cpp:535
T Max() const
Definition: array.cpp:90
void Assemble(int skip_zeros=1)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int Dimension() const
Definition: mesh.hpp:475
void SetTime(double t)
Set physical time (for time-dependent simulations)
void RegisterVisItFields(VisItDataCollection &visit_dc)
void SetMaxIter(int max_iter)
Definition: hypre.cpp:1922
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:133
SurfaceCurrent(ParFiniteElementSpace &H1FESpace, ParFiniteElementSpace &HCurlFESpace, ParDiscreteGradOperator &Grad, Array< int > &kbcs, Array< int > &vbcs, Vector &vbcv)
int SpaceDimension() const
Definition: mesh.hpp:476
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:58
Array< int > bdr_attributes
Definition: mesh.hpp:140
MPI_Comm GetComm() const
Definition: pmesh.hpp:86
PCG solver in hypre.
Definition: hypre.hpp:558
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:323
Class for parallel bilinear form.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:1937
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void ParallelEliminateEssentialBC(const Array< int > &bdr_attr_is_ess, HypreParMatrix &A, const HypreParVector &X, HypreParVector &B) const
void ProjectBdrCoefficientTangent(VectorCoefficient &vcoeff, Array< int > &bdr_attr)
Definition: gridfunc.cpp:1454
Class for parallel bilinear form.
void ComputeSurfaceCurrent(ParGridFunction &k)
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:185
Vector data type.
Definition: vector.hpp:33
Class for parallel grid function.
Definition: pgridfunc.hpp:31
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:143
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:1958
Class for parallel meshes.
Definition: pmesh.hpp:28
void SetSingularProblem()
Set this option when solving a curl-curl problem with zero mass term.
Definition: hypre.hpp:779
Integrator for (Q u, v) for VectorFiniteElements.
Definition: bilininteg.hpp:458
virtual void Mult(const Vector &x, Vector &y) const
Operator application.
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Definition: gridfunc.hpp:170
void ParallelProject(Vector &tv) const
Returns the vector restricted to the true dofs.
Definition: pgridfunc.cpp:120