MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex25p.cpp
Go to the documentation of this file.
1 // MFEM Example 25 - Parallel Version
2 //
3 // Compile with: make ex25p
4 //
5 // Sample runs: mpirun -np 4 ex25p -o 2 -f 1.0 -rs 1 -rp 1 -prob 0
6 // mpirun -np 4 ex25p -o 3 -f 10.0 -rs 1 -rp 1 -prob 1
7 // mpirun -np 4 ex25p -o 3 -f 5.0 -rs 3 -rp 1 -prob 2
8 // mpirun -np 4 ex25p -o 2 -f 1.0 -rs 1 -rp 1 -prob 3
9 // mpirun -np 4 ex25p -o 2 -f 1.0 -rs 2 -rp 2 -prob 0 -m ../data/beam-quad.mesh
10 // mpirun -np 4 ex25p -o 2 -f 8.0 -rs 2 -rp 2 -prob 4 -m ../data/inline-quad.mesh
11 // mpirun -np 4 ex25p -o 2 -f 2.0 -rs 1 -rp 1 -prob 4 -m ../data/inline-hex.mesh
12 //
13 // Device sample runs:
14 // mpirun -np 4 ex25p -o 1 -f 3.0 -rs 3 -rp 1 -prob 2 -pa -d cuda
15 // mpirun -np 4 ex25p -o 2 -f 1.0 -rs 1 -rp 1 -prob 3 -pa -d cuda
16 //
17 // Description: This example code solves a simple electromagnetic wave
18 // propagation problem corresponding to the second order
19 // indefinite Maxwell equation
20 //
21 // (1/mu) * curl curl E - \omega^2 * epsilon E = f
22 //
23 // with a Perfectly Matched Layer (PML).
24 //
25 // The example demonstrates discretization with Nedelec finite
26 // elements in 2D or 3D, as well as the use of complex-valued
27 // bilinear and linear forms. Several test problems are included,
28 // with prob = 0-3 having known exact solutions, see "On perfectly
29 // matched layers for discontinuous Petrov-Galerkin methods" by
30 // Vaziri Astaneh, Keith, Demkowicz, Comput Mech 63, 2019.
31 //
32 // We recommend viewing Example 22 before viewing this example.
33 
34 #include "mfem.hpp"
35 #include <fstream>
36 #include <iostream>
37 
38 #ifdef _WIN32
39 #define jn(n, x) _jn(n, x)
40 #define yn(n, x) _yn(n, x)
41 #endif
42 
43 using namespace std;
44 using namespace mfem;
45 
46 // Class for setting up a simple Cartesian PML region
47 class CartesianPML
48 {
49 private:
50  Mesh *mesh;
51 
52  int dim;
53 
54  // Length of the PML Region in each direction
55  Array2D<double> length;
56 
57  // Computational Domain Boundary
58  Array2D<double> comp_dom_bdr;
59 
60  // Domain Boundary
61  Array2D<double> dom_bdr;
62 
63  // Integer Array identifying elements in the PML
64  // 0: in the PML, 1: not in the PML
65  Array<int> elems;
66 
67  // Compute Domain and Computational Domain Boundaries
68  void SetBoundaries();
69 
70 public:
71  // Constructor
72  CartesianPML(Mesh *mesh_,Array2D<double> length_);
73 
74  // Return Computational Domain Boundary
75  Array2D<double> GetCompDomainBdr() {return comp_dom_bdr;}
76 
77  // Return Domain Boundary
78  Array2D<double> GetDomainBdr() {return dom_bdr;}
79 
80  // Return Markers list for elements
81  Array<int> * GetMarkedPMLElements() {return &elems;}
82 
83  // Mark elements in the PML region
84  void SetAttributes(ParMesh *pmesh);
85 
86  // PML complex stretching function
87  void StretchFunction(const Vector &x, vector<complex<double>> &dxs);
88 };
89 
90 // Class for returning the PML coefficients of the bilinear form
91 class PMLDiagMatrixCoefficient : public VectorCoefficient
92 {
93 private:
94  CartesianPML * pml = nullptr;
95  void (*Function)(const Vector &, CartesianPML *, Vector &);
96 public:
97  PMLDiagMatrixCoefficient(int dim, void(*F)(const Vector &, CartesianPML *,
98  Vector &),
99  CartesianPML * pml_)
100  : VectorCoefficient(dim), pml(pml_), Function(F)
101  {}
102 
103  using VectorCoefficient::Eval;
104 
105  virtual void Eval(Vector &K, ElementTransformation &T,
106  const IntegrationPoint &ip)
107  {
108  double x[3];
109  Vector transip(x, 3);
110  T.Transform(ip, transip);
111  K.SetSize(vdim);
112  (*Function)(transip, pml, K);
113  }
114 };
115 
116 void maxwell_solution(const Vector &x, vector<complex<double>> &Eval);
117 
118 void E_bdr_data_Re(const Vector &x, Vector &E);
119 void E_bdr_data_Im(const Vector &x, Vector &E);
120 
121 void E_exact_Re(const Vector &x, Vector &E);
122 void E_exact_Im(const Vector &x, Vector &E);
123 
124 void source(const Vector &x, Vector & f);
125 
126 // Functions for computing the necessary coefficients after PML stretching.
127 // J is the Jacobian matrix of the stretching function
128 void detJ_JT_J_inv_Re(const Vector &x, CartesianPML * pml, Vector & D);
129 void detJ_JT_J_inv_Im(const Vector &x, CartesianPML * pml, Vector & D);
130 void detJ_JT_J_inv_abs(const Vector &x, CartesianPML * pml, Vector & D);
131 
132 void detJ_inv_JT_J_Re(const Vector &x, CartesianPML * pml, Vector & D);
133 void detJ_inv_JT_J_Im(const Vector &x, CartesianPML * pml, Vector & D);
134 void detJ_inv_JT_J_abs(const Vector &x, CartesianPML * pml, Vector & D);
135 
138 
139 double mu = 1.0;
140 double epsilon = 1.0;
141 double omega;
142 int dim;
143 bool exact_known = false;
144 
146 {
147  beam, // Wave propagating in a beam-like domain
148  disc, // Point source propagating in the square-disc domain
149  lshape, // Point source propagating in the L-shape domain
150  fichera, // Point source propagating in the fichera domain
151  load_src // Approximated point source with PML all around
152 };
154 
155 int main(int argc, char *argv[])
156 {
157  // 1. Initialize MPI and HYPRE.
158  Mpi::Init(argc, argv);
159  int num_procs = Mpi::WorldSize();
160  int myid = Mpi::WorldRank();
161  Hypre::Init();
162 
163  // 2. Parse command-line options.
164  const char *mesh_file = nullptr;
165  int order = 1;
166  int ref_levels = 1;
167  int par_ref_levels = 2;
168  int iprob = 4;
169  double freq = 5.0;
170  bool herm_conv = true;
171  bool slu_solver = false;
172  bool mumps_solver = false;
173  bool visualization = 1;
174  bool pa = false;
175  const char *device_config = "cpu";
176 
177  OptionsParser args(argc, argv);
178  args.AddOption(&mesh_file, "-m", "--mesh",
179  "Mesh file to use.");
180  args.AddOption(&order, "-o", "--order",
181  "Finite element order (polynomial degree).");
182  args.AddOption(&iprob, "-prob", "--problem", "Problem case"
183  " 0: beam, 1: disc, 2: lshape, 3: fichera, 4: General");
184  args.AddOption(&ref_levels, "-rs", "--refinements-serial",
185  "Number of serial refinements");
186  args.AddOption(&par_ref_levels, "-rp", "--refinements-parallel",
187  "Number of parallel refinements");
188  args.AddOption(&mu, "-mu", "--permeability",
189  "Permeability of free space (or 1/(spring constant)).");
190  args.AddOption(&epsilon, "-eps", "--permittivity",
191  "Permittivity of free space (or mass constant).");
192  args.AddOption(&freq, "-f", "--frequency",
193  "Frequency (in Hz).");
194  args.AddOption(&herm_conv, "-herm", "--hermitian", "-no-herm",
195  "--no-hermitian", "Use convention for Hermitian operators.");
196 #ifdef MFEM_USE_SUPERLU
197  args.AddOption(&slu_solver, "-slu", "--superlu", "-no-slu",
198  "--no-superlu", "Use the SuperLU Solver.");
199 #endif
200 #ifdef MFEM_USE_MUMPS
201  args.AddOption(&mumps_solver, "-mumps", "--mumps-solver", "-no-mumps",
202  "--no-mumps-solver", "Use the MUMPS Solver.");
203 #endif
204  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
205  "--no-visualization",
206  "Enable or disable GLVis visualization.");
207  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
208  "--no-partial-assembly", "Enable Partial Assembly.");
209  args.AddOption(&device_config, "-d", "--device",
210  "Device configuration string, see Device::Configure().");
211  args.Parse();
212  if (slu_solver && mumps_solver)
213  {
214  if (myid == 0)
215  cout << "WARNING: Both SuperLU and MUMPS have been selected,"
216  << " please choose either one." << endl
217  << " Defaulting to SuperLU." << endl;
218  mumps_solver = false;
219  }
220 
221  if (iprob > 4) { iprob = 4; }
222  prob = (prob_type)iprob;
223 
224  // 3. Enable hardware devices such as GPUs, and programming models such as
225  // CUDA, OCCA, RAJA and OpenMP based on command line options.
226  Device device(device_config);
227  if (myid == 0) { device.Print(); }
228 
229  // 4. Setup the (serial) mesh on all processors.
230  if (!mesh_file)
231  {
232  exact_known = true;
233  switch (prob)
234  {
235  case beam:
236  mesh_file = "../data/beam-hex.mesh";
237  break;
238  case disc:
239  mesh_file = "../data/square-disc.mesh";
240  break;
241  case lshape:
242  mesh_file = "../data/l-shape.mesh";
243  break;
244  case fichera:
245  mesh_file = "../data/fichera.mesh";
246  break;
247  default:
248  exact_known = false;
249  mesh_file = "../data/inline-quad.mesh";
250  break;
251  }
252  }
253 
254  if (!args.Good())
255  {
256  if (myid == 0)
257  {
258  args.PrintUsage(cout);
259  }
260  return 1;
261  }
262  if (myid == 0)
263  {
264  args.PrintOptions(cout);
265  }
266 
267  Mesh * mesh = new Mesh(mesh_file, 1, 1);
268  dim = mesh->Dimension();
269 
270  // Angular frequency
271  omega = 2.0 * M_PI * freq;
272 
273  // Setup PML length
274  Array2D<double> length(dim, 2); length = 0.0;
275 
276  // 5. Setup the Cartesian PML region.
277  switch (prob)
278  {
279  case disc:
280  length = 0.2;
281  break;
282  case lshape:
283  length(0, 0) = 0.1;
284  length(1, 0) = 0.1;
285  break;
286  case fichera:
287  length(0, 1) = 0.5;
288  length(1, 1) = 0.5;
289  length(2, 1) = 0.5;
290  break;
291  case beam:
292  length(0, 1) = 2.0;
293  break;
294  default:
295  length = 0.25;
296  break;
297  }
298  CartesianPML * pml = new CartesianPML(mesh,length);
299  comp_domain_bdr = pml->GetCompDomainBdr();
300  domain_bdr = pml->GetDomainBdr();
301 
302  // 6. Refine the serial mesh on all processors to increase the resolution.
303  for (int l = 0; l < ref_levels; l++)
304  {
305  mesh->UniformRefinement();
306  }
307 
308  // 7. Define a parallel mesh by a partitioning of the serial mesh.
309  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
310  delete mesh;
311  {
312  for (int l = 0; l < par_ref_levels; l++)
313  {
314  pmesh->UniformRefinement();
315  }
316  }
317 
318  // 8. Set element attributes in order to distinguish elements in the PML
319  pml->SetAttributes(pmesh);
320 
321  // 9. Define a parallel finite element space on the parallel mesh. Here we
322  // use the Nedelec finite elements of the specified order.
323  FiniteElementCollection *fec = new ND_FECollection(order, dim);
324  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
325  HYPRE_BigInt size = fespace->GlobalTrueVSize();
326  if (myid == 0)
327  {
328  cout << "Number of finite element unknowns: " << size << endl;
329  }
330 
331  // 10. Determine the list of true (i.e. parallel conforming) essential
332  // boundary dofs. In this example, the boundary conditions are defined
333  // based on the specific mesh and the problem type.
334  Array<int> ess_tdof_list;
335  Array<int> ess_bdr;
336  if (pmesh->bdr_attributes.Size())
337  {
338  ess_bdr.SetSize(pmesh->bdr_attributes.Max());
339  ess_bdr = 1;
340  if (prob == lshape || prob == fichera)
341  {
342  ess_bdr = 0;
343  for (int j = 0; j < pmesh->GetNBE(); j++)
344  {
345  Vector center(dim);
346  int bdrgeom = pmesh->GetBdrElementBaseGeometry(j);
348  tr->Transform(Geometries.GetCenter(bdrgeom),center);
349  int k = pmesh->GetBdrAttribute(j);
350  switch (prob)
351  {
352  case lshape:
353  if (center[0] == 1.0 || center[0] == 0.5 || center[1] == 0.5)
354  {
355  ess_bdr[k - 1] = 1;
356  }
357  break;
358  case fichera:
359  if (center[0] == -1.0 || center[0] == 0.0 ||
360  center[1] == 0.0 || center[2] == 0.0)
361  {
362  ess_bdr[k - 1] = 1;
363  }
364  break;
365  default:
366  break;
367  }
368  }
369  }
370  }
371  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
372 
373  // 11. Setup Complex Operator convention
375  herm_conv ? ComplexOperator::HERMITIAN : ComplexOperator::BLOCK_SYMMETRIC;
376 
377  // 12. Set up the parallel linear form b(.) which corresponds to the
378  // right-hand side of the FEM linear system.
380  ParComplexLinearForm b(fespace, conv);
381  if (prob == load_src)
382  {
384  }
385  b.Vector::operator=(0.0);
386  b.Assemble();
387 
388  // 13. Define the solution vector x as a parallel complex finite element grid
389  // function corresponding to fespace.
390  ParComplexGridFunction x(fespace);
391  x = 0.0;
394  x.ProjectBdrCoefficientTangent(E_Re, E_Im, ess_bdr);
395 
396  // 14. Set up the parallel sesquilinear form a(.,.)
397  //
398  // In Comp
399  // Domain: 1/mu (Curl E, Curl F) - omega^2 * epsilon (E,F)
400  //
401  // In PML: 1/mu (1/det(J) J^T J Curl E, Curl F)
402  // - omega^2 * epsilon (det(J) * (J^T J)^-1 * E, F)
403  //
404  // where J denotes the Jacobian Matrix of the PML Stretching function
405  Array<int> attr;
406  Array<int> attrPML;
407  if (pmesh->attributes.Size())
408  {
409  attr.SetSize(pmesh->attributes.Max());
410  attrPML.SetSize(pmesh->attributes.Max());
411  attr = 0; attr[0] = 1;
412  attrPML = 0;
413  if (pmesh->attributes.Max() > 1)
414  {
415  attrPML[1] = 1;
416  }
417  }
418 
419  ConstantCoefficient muinv(1.0/mu);
420  ConstantCoefficient omeg(-pow(omega, 2) * epsilon);
421  RestrictedCoefficient restr_muinv(muinv,attr);
422  RestrictedCoefficient restr_omeg(omeg,attr);
423 
424  // Integrators inside the computational domain (excluding the PML region)
425  ParSesquilinearForm a(fespace, conv);
426  a.AddDomainIntegrator(new CurlCurlIntegrator(restr_muinv),NULL);
427  a.AddDomainIntegrator(new VectorFEMassIntegrator(restr_omeg),NULL);
428 
429  int cdim = (dim == 2) ? 1 : dim;
430  PMLDiagMatrixCoefficient pml_c1_Re(cdim,detJ_inv_JT_J_Re, pml);
431  PMLDiagMatrixCoefficient pml_c1_Im(cdim,detJ_inv_JT_J_Im, pml);
432  ScalarVectorProductCoefficient c1_Re(muinv,pml_c1_Re);
433  ScalarVectorProductCoefficient c1_Im(muinv,pml_c1_Im);
434  VectorRestrictedCoefficient restr_c1_Re(c1_Re,attrPML);
435  VectorRestrictedCoefficient restr_c1_Im(c1_Im,attrPML);
436 
437  PMLDiagMatrixCoefficient pml_c2_Re(dim, detJ_JT_J_inv_Re,pml);
438  PMLDiagMatrixCoefficient pml_c2_Im(dim, detJ_JT_J_inv_Im,pml);
439  ScalarVectorProductCoefficient c2_Re(omeg,pml_c2_Re);
440  ScalarVectorProductCoefficient c2_Im(omeg,pml_c2_Im);
441  VectorRestrictedCoefficient restr_c2_Re(c2_Re,attrPML);
442  VectorRestrictedCoefficient restr_c2_Im(c2_Im,attrPML);
443 
444  // Integrators inside the PML region
445  a.AddDomainIntegrator(new CurlCurlIntegrator(restr_c1_Re),
446  new CurlCurlIntegrator(restr_c1_Im));
447  a.AddDomainIntegrator(new VectorFEMassIntegrator(restr_c2_Re),
448  new VectorFEMassIntegrator(restr_c2_Im));
449 
450  // 15. Assemble the parallel bilinear form and the corresponding linear
451  // system, applying any necessary transformations such as: parallel
452  // assembly, eliminating boundary conditions, applying conforming
453  // constraints for non-conforming AMR, etc.
454  if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
455  a.Assemble();
456 
457  OperatorPtr Ah;
458  Vector B, X;
459  a.FormLinearSystem(ess_tdof_list, x, b, Ah, X, B);
460 
461  // 16. Solve using a direct or an iterative solver
462 #ifdef MFEM_USE_SUPERLU
463  if (!pa && slu_solver)
464  {
465  // Transform to monolithic HypreParMatrix
466  HypreParMatrix *A = Ah.As<ComplexHypreParMatrix>()->GetSystemMatrix();
467  SuperLURowLocMatrix SA(*A);
468  SuperLUSolver superlu(MPI_COMM_WORLD);
469  superlu.SetPrintStatistics(false);
470  superlu.SetSymmetricPattern(false);
472  superlu.SetOperator(SA);
473  superlu.Mult(B, X);
474  delete A;
475  }
476 #endif
477 #ifdef MFEM_USE_MUMPS
478  if (!pa && mumps_solver)
479  {
480  HypreParMatrix *A = Ah.As<ComplexHypreParMatrix>()->GetSystemMatrix();
481  MUMPSSolver mumps;
482  mumps.SetPrintLevel(0);
483  mumps.SetMatrixSymType(MUMPSSolver::MatType::UNSYMMETRIC);
484  mumps.SetOperator(*A);
485  mumps.Mult(B,X);
486  delete A;
487  }
488 #endif
489  // 16a. Set up the parallel Bilinear form a(.,.) for the preconditioner
490  //
491  // In Comp
492  // Domain: 1/mu (Curl E, Curl F) + omega^2 * epsilon (E,F)
493  //
494  // In PML: 1/mu (abs(1/det(J) J^T J) Curl E, Curl F)
495  // + omega^2 * epsilon (abs(det(J) * (J^T J)^-1) * E, F)
496  if (pa || (!slu_solver && !mumps_solver))
497  {
498  ConstantCoefficient absomeg(pow(omega, 2) * epsilon);
499  RestrictedCoefficient restr_absomeg(absomeg,attr);
500 
501  ParBilinearForm prec(fespace);
502  prec.AddDomainIntegrator(new CurlCurlIntegrator(restr_muinv));
503  prec.AddDomainIntegrator(new VectorFEMassIntegrator(restr_absomeg));
504 
505  PMLDiagMatrixCoefficient pml_c1_abs(cdim,detJ_inv_JT_J_abs, pml);
506  ScalarVectorProductCoefficient c1_abs(muinv,pml_c1_abs);
507  VectorRestrictedCoefficient restr_c1_abs(c1_abs,attrPML);
508 
509  PMLDiagMatrixCoefficient pml_c2_abs(dim, detJ_JT_J_inv_abs,pml);
510  ScalarVectorProductCoefficient c2_abs(absomeg,pml_c2_abs);
511  VectorRestrictedCoefficient restr_c2_abs(c2_abs,attrPML);
512 
513  prec.AddDomainIntegrator(new CurlCurlIntegrator(restr_c1_abs));
514  prec.AddDomainIntegrator(new VectorFEMassIntegrator(restr_c2_abs));
515 
516  if (pa) { prec.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
517  prec.Assemble();
518 
519  // 16b. Define and apply a parallel GMRES solver for AU=B with a block
520  // diagonal preconditioner based on hypre's AMS preconditioner.
521  Array<int> offsets(3);
522  offsets[0] = 0;
523  offsets[1] = fespace->GetTrueVSize();
524  offsets[2] = fespace->GetTrueVSize();
525  offsets.PartialSum();
526 
527  Operator *pc_r = nullptr;
528  Operator *pc_i = nullptr;
529  int s = (conv == ComplexOperator::HERMITIAN) ? -1.0 : 1.0;
530  if (pa)
531  {
532  // Jacobi Smoother
533  OperatorJacobiSmoother *d00 = new OperatorJacobiSmoother(prec, ess_tdof_list);
534  ScaledOperator *d11 = new ScaledOperator(d00, s);
535  pc_r = d00;
536  pc_i = d11;
537  }
538  else
539  {
540  OperatorPtr PCOpAh;
541  prec.FormSystemMatrix(ess_tdof_list, PCOpAh);
542 
543  // Hypre AMS
544  HypreAMS *ams00 = new HypreAMS(*PCOpAh.As<HypreParMatrix>(), fespace);
545  ScaledOperator *ams11 = new ScaledOperator(ams00, s);
546  pc_r = ams00;
547  pc_i = ams11;
548  }
549 
550  BlockDiagonalPreconditioner BlockDP(offsets);
551  BlockDP.SetDiagonalBlock(0, pc_r);
552  BlockDP.SetDiagonalBlock(1, pc_i);
553 
554  GMRESSolver gmres(MPI_COMM_WORLD);
555  gmres.SetPrintLevel(1);
556  gmres.SetKDim(200);
557  gmres.SetMaxIter(pa ? 5000 : 2000);
558  gmres.SetRelTol(1e-5);
559  gmres.SetAbsTol(0.0);
560  gmres.SetOperator(*Ah);
561  gmres.SetPreconditioner(BlockDP);
562  gmres.Mult(B, X);
563  }
564 
565  // 17. Recover the parallel grid function corresponding to X. This is the
566  // local finite element solution on each processor.
567  a.RecoverFEMSolution(X, b, x);
568 
569  // If exact is known compute the error
570  if (exact_known)
571  {
574  int order_quad = max(2, 2 * order + 1);
575  const IntegrationRule *irs[Geometry::NumGeom];
576  for (int i = 0; i < Geometry::NumGeom; ++i)
577  {
578  irs[i] = &(IntRules.Get(i, order_quad));
579  }
580 
581  double L2Error_Re = x.real().ComputeL2Error(E_ex_Re, irs,
582  pml->GetMarkedPMLElements());
583  double L2Error_Im = x.imag().ComputeL2Error(E_ex_Im, irs,
584  pml->GetMarkedPMLElements());
585 
586  ParComplexGridFunction x_gf0(fespace);
587  x_gf0 = 0.0;
588  double norm_E_Re, norm_E_Im;
589  norm_E_Re = x_gf0.real().ComputeL2Error(E_ex_Re, irs,
590  pml->GetMarkedPMLElements());
591  norm_E_Im = x_gf0.imag().ComputeL2Error(E_ex_Im, irs,
592  pml->GetMarkedPMLElements());
593 
594  if (myid == 0)
595  {
596  cout << "\n Relative Error (Re part): || E_h - E || / ||E|| = "
597  << L2Error_Re / norm_E_Re
598  << "\n Relative Error (Im part): || E_h - E || / ||E|| = "
599  << L2Error_Im / norm_E_Im
600  << "\n Total Error: "
601  << sqrt(L2Error_Re*L2Error_Re + L2Error_Im*L2Error_Im) << "\n\n";
602  }
603  }
604 
605  // 18. Save the refined mesh and the solution in parallel. This output can be
606  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
607  {
608  ostringstream mesh_name, sol_r_name, sol_i_name;
609  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
610  sol_r_name << "ex25p-sol_r." << setfill('0') << setw(6) << myid;
611  sol_i_name << "ex25p-sol_i." << setfill('0') << setw(6) << myid;
612 
613  ofstream mesh_ofs(mesh_name.str().c_str());
614  mesh_ofs.precision(8);
615  pmesh->Print(mesh_ofs);
616 
617  ofstream sol_r_ofs(sol_r_name.str().c_str());
618  ofstream sol_i_ofs(sol_i_name.str().c_str());
619  sol_r_ofs.precision(8);
620  sol_i_ofs.precision(8);
621  x.real().Save(sol_r_ofs);
622  x.imag().Save(sol_i_ofs);
623  }
624 
625  // 19. Send the solution by socket to a GLVis server.
626  if (visualization)
627  {
628  // Define visualization keys for GLVis (see GLVis documentation)
629  string keys;
630  keys = (dim == 3) ? "keys macF\n" : keys = "keys amrRljcUUuu\n";
631  if (prob == beam && dim == 3) {keys = "keys macFFiYYYYYYYYYYYYYYYYYY\n";}
632  if (prob == beam && dim == 2) {keys = "keys amrRljcUUuuu\n"; }
633 
634  char vishost[] = "localhost";
635  int visport = 19916;
636 
637  {
638  socketstream sol_sock_re(vishost, visport);
639  sol_sock_re.precision(8);
640  sol_sock_re << "parallel " << num_procs << " " << myid << "\n"
641  << "solution\n" << *pmesh << x.real() << keys
642  << "window_title 'Solution real part'" << flush;
643  MPI_Barrier(MPI_COMM_WORLD); // try to prevent streams from mixing
644  }
645 
646  {
647  socketstream sol_sock_im(vishost, visport);
648  sol_sock_im.precision(8);
649  sol_sock_im << "parallel " << num_procs << " " << myid << "\n"
650  << "solution\n" << *pmesh << x.imag() << keys
651  << "window_title 'Solution imag part'" << flush;
652  MPI_Barrier(MPI_COMM_WORLD); // try to prevent streams from mixing
653  }
654 
655  {
656  ParGridFunction x_t(fespace);
657  x_t = x.real();
658 
659  socketstream sol_sock(vishost, visport);
660  sol_sock.precision(8);
661  sol_sock << "parallel " << num_procs << " " << myid << "\n"
662  << "solution\n" << *pmesh << x_t << keys << "autoscale off\n"
663  << "window_title 'Harmonic Solution (t = 0.0 T)'"
664  << "pause\n" << flush;
665 
666  if (myid == 0)
667  {
668  cout << "GLVis visualization paused."
669  << " Press space (in the GLVis window) to resume it.\n";
670  }
671 
672  int num_frames = 32;
673  int i = 0;
674  while (sol_sock)
675  {
676  double t = (double)(i % num_frames) / num_frames;
677  ostringstream oss;
678  oss << "Harmonic Solution (t = " << t << " T)";
679 
680  add(cos(2.0*M_PI*t), x.real(), sin(2.0*M_PI*t), x.imag(), x_t);
681  sol_sock << "parallel " << num_procs << " " << myid << "\n";
682  sol_sock << "solution\n" << *pmesh << x_t
683  << "window_title '" << oss.str() << "'" << flush;
684  i++;
685  }
686  }
687  }
688 
689  // 20. Free the used memory.
690  delete pml;
691  delete fespace;
692  delete fec;
693  delete pmesh;
694  return 0;
695 }
696 
697 void source(const Vector &x, Vector &f)
698 {
699  Vector center(dim);
700  double r = 0.0;
701  for (int i = 0; i < dim; ++i)
702  {
703  center(i) = 0.5 * (comp_domain_bdr(i, 0) + comp_domain_bdr(i, 1));
704  r += pow(x[i] - center[i], 2.);
705  }
706  double n = 5.0 * omega * sqrt(epsilon * mu) / M_PI;
707  double coeff = pow(n, 2) / M_PI;
708  double alpha = -pow(n, 2) * r;
709  f = 0.0;
710  f[0] = coeff * exp(alpha);
711 }
712 
713 void maxwell_solution(const Vector &x, vector<complex<double>> &E)
714 {
715  // Initialize
716  for (int i = 0; i < dim; ++i)
717  {
718  E[i] = 0.0;
719  }
720 
721  complex<double> zi = complex<double>(0., 1.);
722  double k = omega * sqrt(epsilon * mu);
723  switch (prob)
724  {
725  case disc:
726  case lshape:
727  case fichera:
728  {
729  Vector shift(dim);
730  shift = 0.0;
731  if (prob == fichera) { shift = 1.0; }
732  if (prob == disc) { shift = -0.5; }
733  if (prob == lshape) { shift = -1.0; }
734 
735  if (dim == 2)
736  {
737  double x0 = x(0) + shift(0);
738  double x1 = x(1) + shift(1);
739  double r = sqrt(x0 * x0 + x1 * x1);
740  double beta = k * r;
741 
742  // Bessel functions
743  complex<double> Ho, Ho_r, Ho_rr;
744  Ho = jn(0, beta) + zi * yn(0, beta);
745  Ho_r = -k * (jn(1, beta) + zi * yn(1, beta));
746  Ho_rr = -k * k * (1.0 / beta *
747  (jn(1, beta) + zi * yn(1, beta)) -
748  (jn(2, beta) + zi * yn(2, beta)));
749 
750  // First derivatives
751  double r_x = x0 / r;
752  double r_y = x1 / r;
753  double r_xy = -(r_x / r) * r_y;
754  double r_xx = (1.0 / r) * (1.0 - r_x * r_x);
755 
756  complex<double> val, val_xx, val_xy;
757  val = 0.25 * zi * Ho;
758  val_xx = 0.25 * zi * (r_xx * Ho_r + r_x * r_x * Ho_rr);
759  val_xy = 0.25 * zi * (r_xy * Ho_r + r_x * r_y * Ho_rr);
760  E[0] = zi / k * (k * k * val + val_xx);
761  E[1] = zi / k * val_xy;
762  }
763  else if (dim == 3)
764  {
765  double x0 = x(0) + shift(0);
766  double x1 = x(1) + shift(1);
767  double x2 = x(2) + shift(2);
768  double r = sqrt(x0 * x0 + x1 * x1 + x2 * x2);
769 
770  double r_x = x0 / r;
771  double r_y = x1 / r;
772  double r_z = x2 / r;
773  double r_xx = (1.0 / r) * (1.0 - r_x * r_x);
774  double r_yx = -(r_y / r) * r_x;
775  double r_zx = -(r_z / r) * r_x;
776 
777  complex<double> val, val_r, val_rr;
778  val = exp(zi * k * r) / r;
779  val_r = val / r * (zi * k * r - 1.0);
780  val_rr = val / (r * r) * (-k * k * r * r
781  - 2.0 * zi * k * r + 2.0);
782 
783  complex<double> val_xx, val_yx, val_zx;
784  val_xx = val_rr * r_x * r_x + val_r * r_xx;
785  val_yx = val_rr * r_x * r_y + val_r * r_yx;
786  val_zx = val_rr * r_x * r_z + val_r * r_zx;
787 
788  complex<double> alpha = zi * k / 4.0 / M_PI / k / k;
789  E[0] = alpha * (k * k * val + val_xx);
790  E[1] = alpha * val_yx;
791  E[2] = alpha * val_zx;
792  }
793  break;
794  }
795  case beam:
796  {
797  // T_10 mode
798  if (dim == 3)
799  {
800  double k10 = sqrt(k * k - M_PI * M_PI);
801  E[1] = -zi * k / M_PI * sin(M_PI*x(2))*exp(zi * k10 * x(0));
802  }
803  else if (dim == 2)
804  {
805  E[1] = -zi * k / M_PI * exp(zi * k * x(0));
806  }
807  break;
808  }
809  default:
810  break;
811  }
812 }
813 
814 void E_exact_Re(const Vector &x, Vector &E)
815 {
816  vector<complex<double>> Eval(E.Size());
817  maxwell_solution(x, Eval);
818  for (int i = 0; i < dim; ++i)
819  {
820  E[i] = Eval[i].real();
821  }
822 }
823 
824 void E_exact_Im(const Vector &x, Vector &E)
825 {
826  vector<complex<double>> Eval(E.Size());
827  maxwell_solution(x, Eval);
828  for (int i = 0; i < dim; ++i)
829  {
830  E[i] = Eval[i].imag();
831  }
832 }
833 
834 void E_bdr_data_Re(const Vector &x, Vector &E)
835 {
836  E = 0.0;
837  bool in_pml = false;
838 
839  for (int i = 0; i < dim; ++i)
840  {
841  // check if in PML
842  if (x(i) - comp_domain_bdr(i, 0) < 0.0 ||
843  x(i) - comp_domain_bdr(i, 1) > 0.0)
844  {
845  in_pml = true;
846  break;
847  }
848  }
849  if (!in_pml)
850  {
851  vector<complex<double>> Eval(E.Size());
852  maxwell_solution(x, Eval);
853  for (int i = 0; i < dim; ++i)
854  {
855  E[i] = Eval[i].real();
856  }
857  }
858 }
859 
860 // Define bdr_data solution
861 void E_bdr_data_Im(const Vector &x, Vector &E)
862 {
863  E = 0.0;
864  bool in_pml = false;
865 
866  for (int i = 0; i < dim; ++i)
867  {
868  // check if in PML
869  if (x(i) - comp_domain_bdr(i, 0) < 0.0 ||
870  x(i) - comp_domain_bdr(i, 1) > 0.0)
871  {
872  in_pml = true;
873  break;
874  }
875  }
876  if (!in_pml)
877  {
878  vector<complex<double>> Eval(E.Size());
879  maxwell_solution(x, Eval);
880  for (int i = 0; i < dim; ++i)
881  {
882  E[i] = Eval[i].imag();
883  }
884  }
885 }
886 
887 void detJ_JT_J_inv_Re(const Vector &x, CartesianPML * pml, Vector & D)
888 {
889  vector<complex<double>> dxs(dim);
890  complex<double> det(1.0, 0.0);
891  pml->StretchFunction(x, dxs);
892 
893  for (int i = 0; i < dim; ++i)
894  {
895  det *= dxs[i];
896  }
897 
898  for (int i = 0; i < dim; ++i)
899  {
900  D(i) = (det / pow(dxs[i], 2)).real();
901  }
902 }
903 
904 void detJ_JT_J_inv_Im(const Vector &x, CartesianPML * pml, Vector & D)
905 {
906  vector<complex<double>> dxs(dim);
907  complex<double> det = 1.0;
908  pml->StretchFunction(x, dxs);
909 
910  for (int i = 0; i < dim; ++i)
911  {
912  det *= dxs[i];
913  }
914 
915  for (int i = 0; i < dim; ++i)
916  {
917  D(i) = (det / pow(dxs[i], 2)).imag();
918  }
919 }
920 
921 void detJ_JT_J_inv_abs(const Vector &x, CartesianPML * pml, Vector & D)
922 {
923  vector<complex<double>> dxs(dim);
924  complex<double> det = 1.0;
925  pml->StretchFunction(x, dxs);
926 
927  for (int i = 0; i < dim; ++i)
928  {
929  det *= dxs[i];
930  }
931 
932  for (int i = 0; i < dim; ++i)
933  {
934  D(i) = abs(det / pow(dxs[i], 2));
935  }
936 }
937 
938 void detJ_inv_JT_J_Re(const Vector &x, CartesianPML * pml, Vector & D)
939 {
940  vector<complex<double>> dxs(dim);
941  complex<double> det(1.0, 0.0);
942  pml->StretchFunction(x, dxs);
943 
944  for (int i = 0; i < dim; ++i)
945  {
946  det *= dxs[i];
947  }
948 
949  // in the 2D case the coefficient is scalar 1/det(J)
950  if (dim == 2)
951  {
952  D = (1.0 / det).real();
953  }
954  else
955  {
956  for (int i = 0; i < dim; ++i)
957  {
958  D(i) = (pow(dxs[i], 2) / det).real();
959  }
960  }
961 }
962 
963 void detJ_inv_JT_J_Im(const Vector &x, CartesianPML * pml, Vector & D)
964 {
965  vector<complex<double>> dxs(dim);
966  complex<double> det = 1.0;
967  pml->StretchFunction(x, dxs);
968 
969  for (int i = 0; i < dim; ++i)
970  {
971  det *= dxs[i];
972  }
973 
974  if (dim == 2)
975  {
976  D = (1.0 / det).imag();
977  }
978  else
979  {
980  for (int i = 0; i < dim; ++i)
981  {
982  D(i) = (pow(dxs[i], 2) / det).imag();
983  }
984  }
985 }
986 
987 void detJ_inv_JT_J_abs(const Vector &x, CartesianPML * pml, Vector & D)
988 {
989  vector<complex<double>> dxs(dim);
990  complex<double> det = 1.0;
991  pml->StretchFunction(x, dxs);
992 
993  for (int i = 0; i < dim; ++i)
994  {
995  det *= dxs[i];
996  }
997 
998  if (dim == 2)
999  {
1000  D = abs(1.0 / det);
1001  }
1002  else
1003  {
1004  for (int i = 0; i < dim; ++i)
1005  {
1006  D(i) = abs(pow(dxs[i], 2) / det);
1007  }
1008  }
1009 }
1010 
1011 CartesianPML::CartesianPML(Mesh *mesh_, Array2D<double> length_)
1012  : mesh(mesh_), length(length_)
1013 {
1014  dim = mesh->Dimension();
1015  SetBoundaries();
1016 }
1017 
1018 void CartesianPML::SetBoundaries()
1019 {
1020  comp_dom_bdr.SetSize(dim, 2);
1021  dom_bdr.SetSize(dim, 2);
1022  Vector pmin, pmax;
1023  mesh->GetBoundingBox(pmin, pmax);
1024  for (int i = 0; i < dim; i++)
1025  {
1026  dom_bdr(i, 0) = pmin(i);
1027  dom_bdr(i, 1) = pmax(i);
1028  comp_dom_bdr(i, 0) = dom_bdr(i, 0) + length(i, 0);
1029  comp_dom_bdr(i, 1) = dom_bdr(i, 1) - length(i, 1);
1030  }
1031 }
1032 
1033 void CartesianPML::SetAttributes(ParMesh *pmesh)
1034 {
1035  // Initialize bdr attributes
1036  for (int i = 0; i < pmesh->GetNBE(); ++i)
1037  {
1038  pmesh->GetBdrElement(i)->SetAttribute(i+1);
1039  }
1040 
1041  int nrelem = pmesh->GetNE();
1042 
1043  // Initialize list with 1
1044  elems.SetSize(nrelem);
1045 
1046  // Loop through the elements and identify which of them are in the PML
1047  for (int i = 0; i < nrelem; ++i)
1048  {
1049  elems[i] = 1;
1050  bool in_pml = false;
1051  Element *el = pmesh->GetElement(i);
1052  Array<int> vertices;
1053 
1054  // Initialize attribute
1055  el->SetAttribute(1);
1056  el->GetVertices(vertices);
1057  int nrvert = vertices.Size();
1058 
1059  // Check if any vertex is in the PML
1060  for (int iv = 0; iv < nrvert; ++iv)
1061  {
1062  int vert_idx = vertices[iv];
1063  double *coords = pmesh->GetVertex(vert_idx);
1064  for (int comp = 0; comp < dim; ++comp)
1065  {
1066  if (coords[comp] > comp_dom_bdr(comp, 1) ||
1067  coords[comp] < comp_dom_bdr(comp, 0))
1068  {
1069  in_pml = true;
1070  break;
1071  }
1072  }
1073  }
1074  if (in_pml)
1075  {
1076  elems[i] = 0;
1077  el->SetAttribute(2);
1078  }
1079  }
1080  pmesh->SetAttributes();
1081 }
1082 
1083 void CartesianPML::StretchFunction(const Vector &x,
1084  vector<complex<double>> &dxs)
1085 {
1086  complex<double> zi = complex<double>(0., 1.);
1087 
1088  double n = 2.0;
1089  double c = 5.0;
1090  double coeff;
1091  double k = omega * sqrt(epsilon * mu);
1092 
1093  // Stretch in each direction independently
1094  for (int i = 0; i < dim; ++i)
1095  {
1096  dxs[i] = 1.0;
1097  if (x(i) >= comp_domain_bdr(i, 1))
1098  {
1099  coeff = n * c / k / pow(length(i, 1), n);
1100  dxs[i] = 1.0 + zi * coeff *
1101  abs(pow(x(i) - comp_domain_bdr(i, 1), n - 1.0));
1102  }
1103  if (x(i) <= comp_domain_bdr(i, 0))
1104  {
1105  coeff = n * c / k / pow(length(i, 0), n);
1106  dxs[i] = 1.0 + zi * coeff *
1107  abs(pow(x(i) - comp_domain_bdr(i, 0), n - 1.0));
1108  }
1109  }
1110 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
Definition: ex25.cpp:148
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:474
int GetBdrAttribute(int i) const
Return the attribute of boundary element i.
Definition: mesh.hpp:1461
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1596
const double * GetVertex(int i) const
Return pointer to vertex i&#39;s coordinates.
Definition: mesh.hpp:1005
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
void SetOperator(const Operator &op)
Set the Operator and perform factorization.
Definition: mumps.cpp:30
Base class for vector Coefficients that optionally depend on time and space.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void GetVertices(Array< int > &v) const =0
Returns element&#39;s vertices.
double epsilon
Definition: ex25.cpp:140
void E_exact_Re(const Vector &x, Vector &E)
Definition: ex25.cpp:737
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
int GetNBE() const
Returns number of boundary elements.
Definition: mesh.hpp:931
Integrator for (curl u, curl v) for Nedelec elements.
void SetPrintLevel(int print_lvl)
Set the error print level for MUMPS.
Definition: mumps.cpp:247
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
void maxwell_solution(const Vector &x, vector< complex< double >> &Eval)
Definition: ex25.cpp:636
void AddDomainIntegrator(LinearFormIntegrator *lfi_real, LinearFormIntegrator *lfi_imag)
Adds new Domain Integrator.
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
void detJ_JT_J_inv_abs(const Vector &x, CartesianPML *pml, Vector &D)
Definition: ex25.cpp:844
ParGridFunction & imag()
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void detJ_inv_JT_J_abs(const Vector &x, CartesianPML *pml, Vector &D)
Definition: ex25.cpp:910
void E_bdr_data_Re(const Vector &x, Vector &E)
Definition: ex25.cpp:757
FDualNumber< tbase > exp(const FDualNumber< tbase > &f)
exp([dual number])
Definition: fdual.hpp:499
const IntegrationPoint & GetCenter(int GeomType)
Return the center of the given Geometry::Type, GeomType.
Definition: geom.hpp:71
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:300
void SetSymmetricPattern(bool sym)
Definition: superlu.cpp:423
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:289
ElementTransformation * GetBdrElementTransformation(int i)
Returns the transformation defining the i-th boundary element.
Definition: mesh.cpp:431
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:179
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void Mult(const Vector &x, Vector &y) const
Solve y = Op^{-1} x.
Definition: mumps.cpp:197
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:580
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:971
void SetMatrixSymType(MatType mtype)
Set the matrix type.
Definition: mumps.cpp:252
Derived vector coefficient that has the value of the parent vector where it is active and is zero oth...
bool exact_known
Definition: ex25.cpp:143
void source(const Vector &x, Vector &f)
Definition: ex25.cpp:620
Geometry Geometries
Definition: fe.cpp:49
void E_exact_Im(const Vector &x, Vector &E)
Definition: ex25.cpp:747
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
Array2D< double > domain_bdr
Definition: ex25.cpp:137
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
constexpr char vishost[]
Derived coefficient that takes the value of the parent coefficient for the active attributes and is z...
Definition: ex25.cpp:149
void E_bdr_data_Im(const Vector &x, Vector &E)
Definition: ex25.cpp:784
void SetColumnPermutation(superlu::ColPerm col_perm)
Definition: superlu.cpp:343
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
constexpr int visport
void detJ_inv_JT_J_Im(const Vector &x, CartesianPML *pml, Vector &D)
Definition: ex25.cpp:886
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
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)
Assemble the local matrix.
void SetKDim(int dim)
Set the number of iteration to perform between restarts, default is 50.
Definition: solvers.hpp:510
void SetPrintStatistics(bool print_stat)
Definition: superlu.cpp:325
const Element * GetElement(int i) const
Definition: mesh.hpp:1033
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
prob_type prob
Definition: ex25.cpp:153
FDualNumber< tbase > cos(const FDualNumber< tbase > &f)
cos([dual number])
Definition: fdual.hpp:471
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
A general vector function coefficient.
double omega
Definition: ex25.cpp:141
void SetAbsTol(double atol)
Definition: solvers.hpp:199
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
Array2D< double > comp_domain_bdr
Definition: ex25.cpp:136
void SetRelTol(double rtol)
Definition: solvers.hpp:198
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Definition: fdual.hpp:572
Scaled Operator B: x -&gt; a A(x).
Definition: operator.hpp:695
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
Vector coefficient defined as a product of scalar and vector coefficients.
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
HYPRE_Int HYPRE_BigInt
GMRES method.
Definition: solvers.hpp:497
void SetAttributes() override
Definition: pmesh.cpp:1579
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
void detJ_JT_J_inv_Re(const Vector &x, CartesianPML *pml, Vector &D)
Definition: ex25.cpp:810
double a
Definition: lissajous.cpp:41
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Class for integration point with weight.
Definition: intrules.hpp:25
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:281
double mu
Definition: ex25.cpp:139
FDualNumber< tbase > sqrt(const FDualNumber< tbase > &f)
sqrt([dual number])
Definition: fdual.hpp:600
int dim
Definition: ex24.cpp:53
void detJ_JT_J_inv_Im(const Vector &x, CartesianPML *pml, Vector &D)
Definition: ex25.cpp:827
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
double freq
Definition: ex24.cpp:54
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
void SetAttribute(const int attr)
Set element&#39;s attribute.
Definition: element.hpp:58
Class for parallel bilinear form.
Geometry::Type GetBdrElementBaseGeometry(int i) const
Definition: mesh.hpp:1065
RefCoord t[3]
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
const double alpha
Definition: ex15.cpp:369
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:411
Vector data type.
Definition: vector.hpp:60
virtual void Transform(const IntegrationPoint &, Vector &)=0
Transform integration point from reference coordinates to physical coordinates and store them in the ...
Definition: ex25.cpp:147
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
ParGridFunction & real()
RefCoord s[3]
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4684
Class for parallel grid function.
Definition: pgridfunc.hpp:32
prob_type
Definition: ex25.cpp:145
void detJ_inv_JT_J_Re(const Vector &x, CartesianPML *pml, Vector &D)
Definition: ex25.cpp:861
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
Abstract data type element.
Definition: element.hpp:28
int main()
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:268
const Element * GetBdrElement(int i) const
Definition: mesh.hpp:1037
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
MUMPS: A Parallel Sparse Direct Solver.
Definition: mumps.hpp:34
double f(const Vector &p)
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150