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