MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
volta_solver.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2020, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "volta_solver.hpp"
13 
14 #ifdef MFEM_USE_MPI
15 
16 using namespace std;
17 namespace mfem
18 {
19 
20 using namespace common;
21 
22 namespace electromagnetics
23 {
24 
25 VoltaSolver::VoltaSolver(ParMesh & pmesh, int order,
26  Array<int> & dbcs, Vector & dbcv,
27  Array<int> & nbcs, Vector & nbcv,
28  Coefficient & epsCoef,
29  double (*phi_bc )(const Vector&),
30  double (*rho_src)(const Vector&),
31  void (*p_src )(const Vector&, Vector&),
32  Vector & point_charges)
33  : myid_(0),
34  num_procs_(1),
35  order_(order),
36  pmesh_(&pmesh),
37  dbcs_(&dbcs),
38  dbcv_(&dbcv),
39  nbcs_(&nbcs),
40  nbcv_(&nbcv),
41  visit_dc_(NULL),
42  H1FESpace_(NULL),
43  HCurlFESpace_(NULL),
44  HDivFESpace_(NULL),
45  L2FESpace_(NULL),
46  divEpsGrad_(NULL),
47  h1Mass_(NULL),
48  h1SurfMass_(NULL),
49  hDivMass_(NULL),
50  hCurlHDivEps_(NULL),
51  hCurlHDiv_(NULL),
52  weakDiv_(NULL),
53  rhod_(NULL),
54  l2_vol_int_(NULL),
55  rt_surf_int_(NULL),
56  grad_(NULL),
57  phi_(NULL),
58  rho_src_(NULL),
59  rho_(NULL),
60  sigma_src_(NULL),
61  e_(NULL),
62  d_(NULL),
63  p_src_(NULL),
64  oneCoef_(1.0),
65  epsCoef_(&epsCoef),
66  phiBCCoef_(NULL),
67  rhoCoef_(NULL),
68  pCoef_(NULL),
69  phi_bc_func_(phi_bc),
70  rho_src_func_(rho_src),
71  p_src_func_(p_src),
72  point_charge_params_(point_charges),
73  point_charges_(0)
74 {
75  // Initialize MPI variables
76  MPI_Comm_size(pmesh_->GetComm(), &num_procs_);
77  MPI_Comm_rank(pmesh_->GetComm(), &myid_);
78 
79  // Define compatible parallel finite element spaces on the parallel
80  // mesh. Here we use arbitrary order H1, Nedelec, and Raviart-Thomas finite
81  // elements.
82  H1FESpace_ = new H1_ParFESpace(pmesh_,order,pmesh_->Dimension());
83  HCurlFESpace_ = new ND_ParFESpace(pmesh_,order,pmesh_->Dimension());
84  HDivFESpace_ = new RT_ParFESpace(pmesh_,order,pmesh_->Dimension());
85  L2FESpace_ = new L2_ParFESpace(pmesh_,order-1,pmesh_->Dimension());
86 
87  // Select surface attributes for Dirichlet BCs
88  ess_bdr_.SetSize(pmesh.bdr_attributes.Max());
89  ess_bdr_ = 0; // Deselect all outer surfaces
90  for (int i=0; i<dbcs_->Size(); i++)
91  {
92  ess_bdr_[(*dbcs_)[i]-1] = 1;
93  }
94 
95  // Setup various coefficients
96 
97  // Potential on outer surface
98  if ( phi_bc_func_ != NULL )
99  {
100  phiBCCoef_ = new FunctionCoefficient(*phi_bc_func_);
101  }
102 
103  // Volume Charge Density
104  if ( rho_src_func_ != NULL )
105  {
106  rhoCoef_ = new FunctionCoefficient(rho_src_func_);
107  }
108 
109  // Polarization
110  if ( p_src_func_ != NULL )
111  {
112  pCoef_ = new VectorFunctionCoefficient(pmesh_->SpaceDimension(),
113  p_src_func_);
114  }
115 
116  // Bilinear Forms
117  divEpsGrad_ = new ParBilinearForm(H1FESpace_);
118  divEpsGrad_->AddDomainIntegrator(new DiffusionIntegrator(*epsCoef_));
119 
120  hDivMass_ = new ParBilinearForm(HDivFESpace_);
122 
123  hCurlHDivEps_ = new ParMixedBilinearForm(HCurlFESpace_,HDivFESpace_);
124  hCurlHDivEps_->AddDomainIntegrator(new VectorFEMassIntegrator(*epsCoef_));
125 
126  rhod_ = new ParLinearForm(H1FESpace_);
127 
128  l2_vol_int_ = new ParLinearForm(L2FESpace_);
129  l2_vol_int_->AddDomainIntegrator(new DomainLFIntegrator(oneCoef_));
130 
131  rt_surf_int_ = new ParLinearForm(HDivFESpace_);
133 
134  // Discrete derivative operator
135  grad_ = new ParDiscreteGradOperator(H1FESpace_, HCurlFESpace_);
136  div_ = new ParDiscreteDivOperator(HDivFESpace_, L2FESpace_);
137 
138  // Build grid functions
139  phi_ = new ParGridFunction(H1FESpace_);
140  d_ = new ParGridFunction(HDivFESpace_);
141  e_ = new ParGridFunction(HCurlFESpace_);
142  rho_ = new ParGridFunction(L2FESpace_);
143 
144  if ( point_charge_params_.Size() > 0 )
145  {
146  int dim = pmesh_->Dimension();
147  int npts = point_charge_params_.Size() / (dim + 1);
148  point_charges_.resize(npts);
149 
150  Vector cent(dim);
151  for (int i=0; i<npts; i++)
152  {
153  for (int d=0; d<dim; d++)
154  {
155  cent[d] = point_charge_params_[(dim + 1) * i + d];
156  }
157  double s = point_charge_params_[(dim + 1) * i + dim];
158 
159  point_charges_[i] = new DeltaCoefficient();
160  point_charges_[i]->SetScale(s);
161  point_charges_[i]->SetDeltaCenter(cent);
162 
163  rhod_->AddDomainIntegrator(new DomainLFIntegrator(*point_charges_[i]));
164  }
165  }
166 
167  if ( rho_src_func_ )
168  {
169  rho_src_ = new ParGridFunction(H1FESpace_);
170 
171  h1Mass_ = new ParBilinearForm(H1FESpace_);
172  h1Mass_->AddDomainIntegrator(new MassIntegrator);
173  }
174 
175  if ( p_src_func_ )
176  {
177  p_src_ = new ParGridFunction(HCurlFESpace_);
178 
179  hCurlHDiv_ = new ParMixedBilinearForm(HCurlFESpace_, HDivFESpace_);
181 
182  weakDiv_ = new ParMixedBilinearForm(HCurlFESpace_, H1FESpace_);
184  }
185 
186  if ( nbcs_->Size() > 0 )
187  {
188  sigma_src_ = new ParGridFunction(H1FESpace_);
189 
190  h1SurfMass_ = new ParBilinearForm(H1FESpace_);
191  h1SurfMass_->AddBoundaryIntegrator(new MassIntegrator);
192  }
193 }
194 
196 {
197  delete phiBCCoef_;
198  delete rhoCoef_;
199  delete pCoef_;
200 
201  delete phi_;
202  delete rho_src_;
203  delete rho_;
204  delete rhod_;
205  delete l2_vol_int_;
206  delete rt_surf_int_;
207  delete sigma_src_;
208  delete d_;
209  delete e_;
210  delete p_src_;
211 
212  delete grad_;
213  delete div_;
214 
215  delete divEpsGrad_;
216  delete h1Mass_;
217  delete h1SurfMass_;
218  delete hDivMass_;
219  delete hCurlHDivEps_;
220  delete hCurlHDiv_;
221  delete weakDiv_;
222 
223  delete H1FESpace_;
224  delete HCurlFESpace_;
225  delete HDivFESpace_;
226  delete L2FESpace_;
227 
228  for (unsigned int i=0; i<point_charges_.size(); i++)
229  {
230  delete point_charges_[i];
231  }
232 
233  map<string,socketstream*>::iterator mit;
234  for (mit=socks_.begin(); mit!=socks_.end(); mit++)
235  {
236  delete mit->second;
237  }
238 }
239 
240 HYPRE_Int
242 {
243  return H1FESpace_->GlobalTrueVSize();
244 }
245 
246 void
248 {
249  HYPRE_Int size_h1 = H1FESpace_->GlobalTrueVSize();
250  HYPRE_Int size_nd = HCurlFESpace_->GlobalTrueVSize();
251  HYPRE_Int size_rt = HDivFESpace_->GlobalTrueVSize();
252  HYPRE_Int size_l2 = L2FESpace_->GlobalTrueVSize();
253  if (myid_ == 0)
254  {
255  cout << "Number of H1 unknowns: " << size_h1 << endl;
256  cout << "Number of H(Curl) unknowns: " << size_nd << endl;
257  cout << "Number of H(Div) unknowns: " << size_rt << endl;
258  cout << "Number of L2 unknowns: " << size_l2 << endl;
259  }
260 }
261 
263 {
264  if (myid_ == 0) { cout << "Assembling ... " << flush; }
265 
266  divEpsGrad_->Assemble();
267  divEpsGrad_->Finalize();
268 
269  hDivMass_->Assemble();
270  hDivMass_->Finalize();
271 
272  hCurlHDivEps_->Assemble();
273  hCurlHDivEps_->Finalize();
274 
275  *rhod_ = 0.0;
276  rhod_->Assemble();
277 
278  l2_vol_int_->Assemble();
279  rt_surf_int_->Assemble();
280 
281  grad_->Assemble();
282  grad_->Finalize();
283 
284  div_->Assemble();
285  div_->Finalize();
286 
287  if ( h1Mass_ )
288  {
289  h1Mass_->Assemble();
290  h1Mass_->Finalize();
291  }
292  if ( h1SurfMass_ )
293  {
294  h1SurfMass_->Assemble();
295  h1SurfMass_->Finalize();
296  }
297  if ( hCurlHDiv_ )
298  {
299  hCurlHDiv_->Assemble();
300  hCurlHDiv_->Finalize();
301  }
302  if ( weakDiv_ )
303  {
304  weakDiv_->Assemble();
305  weakDiv_->Finalize();
306  }
307 
308  if (myid_ == 0) { cout << "done." << endl << flush; }
309 }
310 
311 void
313 {
314  if (myid_ == 0) { cout << "Updating ..." << endl; }
315 
316  // Inform the spaces that the mesh has changed
317  // Note: we don't need to interpolate any GridFunctions on the new mesh
318  // so we pass 'false' to skip creation of any transformation matrices.
319  H1FESpace_->Update(false);
320  HCurlFESpace_->Update(false);
321  HDivFESpace_->Update(false);
322  L2FESpace_->Update(false);
323 
324  // Inform the grid functions that the space has changed.
325  phi_->Update();
326  rhod_->Update();
327  l2_vol_int_->Update();
328  rt_surf_int_->Update();
329  d_->Update();
330  e_->Update();
331  rho_->Update();
332  if ( rho_src_ ) { rho_src_->Update(); }
333  if ( sigma_src_ ) { sigma_src_->Update(); }
334  if ( p_src_ ) { p_src_->Update(); }
335 
336  // Inform the bilinear forms that the space has changed.
337  divEpsGrad_->Update();
338  hDivMass_->Update();
339  hCurlHDivEps_->Update();
340 
341  if ( h1Mass_ ) { h1Mass_->Update(); }
342  if ( h1SurfMass_ ) { h1SurfMass_->Update(); }
343  if ( hCurlHDiv_ ) { hCurlHDiv_->Update(); }
344  if ( weakDiv_ ) { weakDiv_->Update(); }
345 
346  // Inform the other objects that the space has changed.
347  grad_->Update();
348  div_->Update();
349 }
350 
351 void
353 {
354  if (myid_ == 0) { cout << "Running solver ... " << endl; }
355 
356  // Initialize the electric potential with its boundary conditions
357  *phi_ = 0.0;
358 
359  if ( dbcs_->Size() > 0 )
360  {
361  if ( phiBCCoef_ )
362  {
363  // Apply gradient boundary condition
364  phi_->ProjectBdrCoefficient(*phiBCCoef_, ess_bdr_);
365  }
366  else
367  {
368  // Apply piecewise constant boundary condition
369  Array<int> dbc_bdr_attr(pmesh_->bdr_attributes.Max());
370  for (int i=0; i<dbcs_->Size(); i++)
371  {
372  ConstantCoefficient voltage((*dbcv_)[i]);
373  dbc_bdr_attr = 0;
374  if ((*dbcs_)[i] <= dbc_bdr_attr.Size())
375  {
376  dbc_bdr_attr[(*dbcs_)[i]-1] = 1;
377  }
378  phi_->ProjectBdrCoefficient(voltage, dbc_bdr_attr);
379  }
380  }
381  }
382 
383  // Initialize the volumetric charge density
384  if ( rho_src_ )
385  {
386  rho_src_->ProjectCoefficient(*rhoCoef_);
387  h1Mass_->AddMult(*rho_src_, *rhod_);
388  }
389 
390  // Initialize the Polarization
391  if ( p_src_ )
392  {
393  p_src_->ProjectCoefficient(*pCoef_);
394  weakDiv_->AddMult(*p_src_, *rhod_);
395  }
396 
397  // Initialize the surface charge density
398  if ( sigma_src_ )
399  {
400  *sigma_src_ = 0.0;
401 
402  Array<int> nbc_bdr_attr(pmesh_->bdr_attributes.Max());
403  for (int i=0; i<nbcs_->Size(); i++)
404  {
405  ConstantCoefficient sigma_coef((*nbcv_)[i]);
406  nbc_bdr_attr = 0;
407  if ((*nbcs_)[i] <= nbc_bdr_attr.Size())
408  {
409  nbc_bdr_attr[(*nbcs_)[i]-1] = 1;
410  }
411  sigma_src_->ProjectBdrCoefficient(sigma_coef, nbc_bdr_attr);
412  }
413  h1SurfMass_->AddMult(*sigma_src_, *rhod_);
414  }
415 
416  // Determine the essential BC degrees of freedom
417  if ( dbcs_->Size() > 0 )
418  {
419  // From user supplied boundary attributes
420  H1FESpace_->GetEssentialTrueDofs(ess_bdr_, ess_bdr_tdofs_);
421  }
422  else
423  {
424  // Use the first DoF on processor zero by default
425  if ( myid_ == 0 )
426  {
427  ess_bdr_tdofs_.SetSize(1);
428  ess_bdr_tdofs_[0] = 0;
429  }
430  }
431 
432  // Apply essential BC and form linear system
433  HypreParMatrix DivEpsGrad;
434  HypreParVector Phi(H1FESpace_);
435  HypreParVector RHS(H1FESpace_);
436 
437  divEpsGrad_->FormLinearSystem(ess_bdr_tdofs_, *phi_, *rhod_, DivEpsGrad,
438  Phi, RHS);
439 
440  // Define and apply a parallel PCG solver for AX=B with the AMG
441  // preconditioner from hypre.
442  HypreBoomerAMG amg(DivEpsGrad);
443  HyprePCG pcg(DivEpsGrad);
444  pcg.SetTol(1e-12);
445  pcg.SetMaxIter(500);
446  pcg.SetPrintLevel(2);
447  pcg.SetPreconditioner(amg);
448  pcg.Mult(RHS, Phi);
449 
450  // Extract the parallel grid function corresponding to the finite
451  // element approximation Phi. This is the local solution on each
452  // processor.
453  divEpsGrad_->RecoverFEMSolution(Phi, *rhod_, *phi_);
454 
455  // Compute the negative Gradient of the solution vector. This is
456  // the magnetic field corresponding to the scalar potential
457  // represented by phi.
458  grad_->Mult(*phi_, *e_); *e_ *= -1.0;
459 
460  // Compute electric displacement (D) from E and P (if present)
461  if (myid_ == 0) { cout << "Computing D ..." << flush; }
462 
463  ParGridFunction ed(HDivFESpace_);
464  hCurlHDivEps_->Mult(*e_, ed);
465  if ( p_src_ )
466  {
467  hCurlHDiv_->AddMult(*p_src_, ed, -1.0);
468  }
469 
470  HypreParMatrix MassHDiv;
471  Vector ED, D;
472 
473  Array<int> dbc_dofs_d;
474  hDivMass_->FormLinearSystem(dbc_dofs_d, *d_, ed, MassHDiv, D, ED);
475 
476  HyprePCG pcgM(MassHDiv);
477  pcgM.SetTol(1e-12);
478  pcgM.SetMaxIter(500);
479  pcgM.SetPrintLevel(0);
480  HypreDiagScale diagM;
481  pcgM.SetPreconditioner(diagM);
482  pcgM.Mult(ED, D);
483 
484  hDivMass_->RecoverFEMSolution(D, ed, *d_);
485 
486  // Compute charge density from rho = Div(D)
487  div_->Mult(*d_, *rho_);
488 
489  if (myid_ == 0) { cout << "done." << flush; }
490 
491  {
492  // Compute total charge as volume integral of rho
493  double charge_rho = (*l2_vol_int_)(*rho_);
494 
495  // Compute total charge as surface integral of D
496  double charge_D = (*rt_surf_int_)(*d_);
497 
498  if (myid_ == 0)
499  {
500  cout << endl << "Total charge: \n"
501  << " Volume integral of charge density: " << charge_rho
502  << "\n Surface integral of dielectric flux: " << charge_D
503  << endl << flush;
504  }
505  }
506 
507  if (myid_ == 0) { cout << "Solver done. " << endl; }
508 }
509 
510 void
512 {
513  if (myid_ == 0) { cout << "Estimating Error ... " << flush; }
514 
515  // Space for the discontinuous (original) flux
516  DiffusionIntegrator flux_integrator(*epsCoef_);
517  L2_FECollection flux_fec(order_, pmesh_->Dimension());
518  // ND_FECollection flux_fec(order_, pmesh_->Dimension());
519  ParFiniteElementSpace flux_fes(pmesh_, &flux_fec, pmesh_->SpaceDimension());
520 
521  // Space for the smoothed (conforming) flux
522  double norm_p = 1;
523  RT_FECollection smooth_flux_fec(order_-1, pmesh_->Dimension());
524  ParFiniteElementSpace smooth_flux_fes(pmesh_, &smooth_flux_fec);
525 
526  L2ZZErrorEstimator(flux_integrator, *phi_,
527  smooth_flux_fes, flux_fes, errors, norm_p);
528 
529  if (myid_ == 0) { cout << "done." << endl; }
530 }
531 
532 void
534 {
535  visit_dc_ = &visit_dc;
536 
537  visit_dc.RegisterField("Phi", phi_);
538  visit_dc.RegisterField("D", d_);
539  visit_dc.RegisterField("E", e_);
540  visit_dc.RegisterField("Rho", rho_);
541  if ( rho_src_ ) { visit_dc.RegisterField("Rho Source", rho_src_); }
542  if ( p_src_ ) { visit_dc.RegisterField("P Source", p_src_); }
543  if ( sigma_src_ ) { visit_dc.RegisterField("Sigma Source", sigma_src_); }
544 }
545 
546 void
548 {
549  if ( visit_dc_ )
550  {
551  if (myid_ == 0) { cout << "Writing VisIt files ..." << flush; }
552 
553  HYPRE_Int prob_size = this->GetProblemSize();
554  visit_dc_->SetCycle(it);
555  visit_dc_->SetTime(prob_size);
556  visit_dc_->Save();
557 
558  if (myid_ == 0) { cout << " done." << endl; }
559  }
560 }
561 
562 void
564 {
565  if ( myid_ == 0 ) { cout << "Opening GLVis sockets." << endl; }
566 
567  socks_["Phi"] = new socketstream;
568  socks_["Phi"]->precision(8);
569 
570  socks_["D"] = new socketstream;
571  socks_["D"]->precision(8);
572 
573  socks_["E"] = new socketstream;
574  socks_["E"]->precision(8);
575 
576  socks_["Rho"] = new socketstream;
577  socks_["Rho"]->precision(8);
578 
579  if ( rho_src_ )
580  {
581  socks_["RhoSrc"] = new socketstream;
582  socks_["RhoSrc"]->precision(8);
583  }
584  if ( p_src_ )
585  {
586  socks_["PSrc"] = new socketstream;
587  socks_["PSrc"]->precision(8);
588  }
589  if ( sigma_src_ )
590  {
591  socks_["SigmaSrc"] = new socketstream;
592  socks_["SigmaSrc"]->precision(8);
593  }
594 }
595 
596 void
598 {
599  if (myid_ == 0) { cout << "Sending data to GLVis ..." << flush; }
600 
601  char vishost[] = "localhost";
602  int visport = 19916;
603 
604  int Wx = 0, Wy = 0; // window position
605  int Ww = 350, Wh = 350; // window size
606  int offx = Ww+10, offy = Wh+45; // window offsets
607 
608  VisualizeField(*socks_["Phi"], vishost, visport,
609  *phi_, "Electric Potential (Phi)", Wx, Wy, Ww, Wh);
610  Wx += offx;
611 
612  VisualizeField(*socks_["E"], vishost, visport,
613  *e_, "Electric Field (E)", Wx, Wy, Ww, Wh);
614  Wx += offx;
615 
616  VisualizeField(*socks_["D"], vishost, visport,
617  *d_, "Electric Displacement (D)", Wx, Wy, Ww, Wh);
618  Wx += offx;
619 
620  VisualizeField(*socks_["Rho"], vishost, visport,
621  *rho_, "Charge Density", Wx, Wy, Ww, Wh);
622  Wx = 0; Wy += offy; // next line
623 
624  if ( rho_src_ )
625  {
626  VisualizeField(*socks_["RhoSrc"], vishost, visport,
627  *rho_src_, "Charge Density Source (Rho)", Wx, Wy, Ww, Wh);
628  Wx += offx;
629  }
630  if ( p_src_ )
631  {
632  VisualizeField(*socks_["PSrc"], vishost, visport,
633  *p_src_, "Electric Polarization Source (P)",
634  Wx, Wy, Ww, Wh);
635  Wx += offx;
636  }
637  if ( sigma_src_ )
638  {
639  VisualizeField(*socks_["SigmaSrc"], vishost, visport,
640  *sigma_src_, "Surface Charge Density Source (Sigma)",
641  Wx, Wy, Ww, Wh);
642  // Wx += offx; // not used
643  }
644  if (myid_ == 0) { cout << " done." << endl; }
645 }
646 
647 } // namespace electromagnetics
648 
649 } // namespace mfem
650 
651 #endif // MFEM_USE_MPI
void SetTol(double tol)
Definition: hypre.cpp:2305
int Size() const
Logical size of the array.
Definition: array.hpp:124
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:770
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
Subclass constant coefficient.
Definition: coefficient.hpp:67
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void AddMult(const Vector &x, Vector &y, const double a=1.0) const
int offx
virtual void Update(bool want_transform=true)
Definition: pfespace.cpp:2868
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
int Wx
Delta function coefficient.
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
virtual void Save()
Save the collection and a VisIt root file.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void RegisterVisItFields(VisItDataCollection &visit_dc)
virtual void ProjectCoefficient(Coefficient &coeff)
Definition: pgridfunc.cpp:301
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2320
void Update(ParFiniteElementSpace *pf=NULL)
Update the object according to the given new FE space *pf.
Definition: plinearform.cpp:21
virtual void Update(FiniteElementSpace *nfes=NULL)
void GetErrorEstimates(Vector &errors)
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
Class for parallel linear form.
Definition: plinearform.hpp:26
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
Definition: pgridfunc.cpp:81
Jacobi preconditioner in hypre.
Definition: hypre.hpp:859
Data collection with VisIt I/O routines.
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
int Dimension() const
Definition: mesh.hpp:744
void SetTime(double t)
Set physical time (for time-dependent simulations)
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2310
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:191
int SpaceDimension() const
Definition: mesh.hpp:745
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:73
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
MPI_Comm GetComm() const
Definition: pmesh.hpp:230
PCG solver in hypre.
Definition: hypre.hpp:743
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:409
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
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)
Definition: pgridfunc.cpp:706
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Class for parallel bilinear form using different test and trial FE spaces.
int Wh
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
void AddBoundaryIntegrator(BilinearFormIntegrator *bfi)
Adds new Boundary Integrator. Assumes ownership of bfi.
int dim
Definition: ex24.cpp:43
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2325
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Class for parallel bilinear form.
class for C-function coefficient
Vector data type.
Definition: vector.hpp:48
void VisualizeField(socketstream &sock, const char *vishost, int visport, GridFunction &gf, const char *title, int x, int y, int w, int h, const char *keys, bool vec)
Definition: fem_extras.cpp:94
int offy
int Wy
int Ww
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2348
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:143