MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex22p.cpp
Go to the documentation of this file.
1 // MFEM Example 22 - Parallel Version
2 //
3 // Compile with: make ex22p
4 //
5 // Sample runs: mpirun -np 4 ex22p -m ../data/inline-segment.mesh -o 3
6 // mpirun -np 4 ex22p -m ../data/inline-tri.mesh -o 3
7 // mpirun -np 4 ex22p -m ../data/inline-quad.mesh -o 3
8 // mpirun -np 4 ex22p -m ../data/inline-quad.mesh -o 3 -p 1
9 // mpirun -np 4 ex22p -m ../data/inline-quad.mesh -o 3 -p 2
10 // mpirun -np 4 ex22p -m ../data/inline-tet.mesh -o 2
11 // mpirun -np 4 ex22p -m ../data/inline-hex.mesh -o 2
12 // mpirun -np 4 ex22p -m ../data/inline-hex.mesh -o 2 -p 1
13 // mpirun -np 4 ex22p -m ../data/inline-hex.mesh -o 2 -p 2
14 // mpirun -np 4 ex22p -m ../data/star.mesh -o 2 -sigma 10.0
15 //
16 // Description: This example code demonstrates the use of MFEM to define and
17 // solve simple complex-valued linear systems. It implements three
18 // variants of a damped harmonic oscillator:
19 //
20 // 1) A scalar H1 field
21 // -Div(a Grad u) - omega^2 b u + i omega c u = 0
22 //
23 // 2) A vector H(Curl) field
24 // Curl(a Curl u) - omega^2 b u + i omega c u = 0
25 //
26 // 3) A vector H(Div) field
27 // -Grad(a Div u) - omega^2 b u + i omega c u = 0
28 //
29 // In each case the field is driven by a forced oscillation, with
30 // angular frequency omega, imposed at the boundary or a portion
31 // of the boundary.
32 //
33 // In electromagnetics the coefficients are typically named the
34 // permeability, mu = 1/a, permittivity, epsilon = b, and
35 // conductivity, sigma = c. The user can specify these constants
36 // using either set of names.
37 //
38 // The example also demonstrates how to display a time-varying
39 // solution as a sequence of fields sent to a single GLVis socket.
40 //
41 // We recommend viewing examples 1, 3 and 4 before viewing this
42 // example.
43 
44 
45 #include "mfem.hpp"
46 #include <fstream>
47 #include <iostream>
48 
49 using namespace std;
50 using namespace mfem;
51 
52 static double mu_ = 1.0;
53 static double epsilon_ = 1.0;
54 static double sigma_ = 20.0;
55 static double omega_ = 10.0;
56 
57 double u0_real_exact(const Vector &);
58 double u0_imag_exact(const Vector &);
59 
60 void u1_real_exact(const Vector &, Vector &);
61 void u1_imag_exact(const Vector &, Vector &);
62 
63 void u2_real_exact(const Vector &, Vector &);
64 void u2_imag_exact(const Vector &, Vector &);
65 
66 bool check_for_inline_mesh(const char * mesh_file);
67 
68 int main(int argc, char *argv[])
69 {
70  // 1. Initialize MPI.
71  int num_procs, myid;
72  MPI_Init(&argc, &argv);
73  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
74  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
75 
76  // 2. Parse command-line options.
77  const char *mesh_file = "../data/inline-quad.mesh";
78  int ser_ref_levels = 1;
79  int par_ref_levels = 1;
80  int order = 1;
81  int prob = 0;
82  double freq = -1.0;
83  double a_coef = 0.0;
84  bool visualization = 1;
85  bool herm_conv = true;
86  bool exact_sol = true;
87 
88  OptionsParser args(argc, argv);
89  args.AddOption(&mesh_file, "-m", "--mesh",
90  "Mesh file to use.");
91  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
92  "Number of times to refine the mesh uniformly in serial.");
93  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
94  "Number of times to refine the mesh uniformly in parallel.");
95  args.AddOption(&order, "-o", "--order",
96  "Finite element order (polynomial degree).");
97  args.AddOption(&prob, "-p", "--problem-type",
98  "Choose between 0: H_1, 1: H(Curl), or 2: H(Div) "
99  "damped harmonic oscillator.");
100  args.AddOption(&a_coef, "-a", "--stiffness-coef",
101  "Stiffness coefficient (spring constant or 1/mu).");
102  args.AddOption(&epsilon_, "-b", "--mass-coef",
103  "Mass coefficient (or epsilon).");
104  args.AddOption(&sigma_, "-c", "--damping-coef",
105  "Damping coefficient (or sigma).");
106  args.AddOption(&mu_, "-mu", "--permeability",
107  "Permeability of free space (or 1/(spring constant)).");
108  args.AddOption(&epsilon_, "-eps", "--permittivity",
109  "Permittivity of free space (or mass constant).");
110  args.AddOption(&sigma_, "-sigma", "--conductivity",
111  "Conductivity (or damping constant).");
112  args.AddOption(&freq, "-f", "--frequency",
113  "Frequency (in Hz).");
114  args.AddOption(&herm_conv, "-herm", "--hermitian", "-no-herm",
115  "--no-hermitian", "Use convention for Hermitian operators.");
116  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
117  "--no-visualization",
118  "Enable or disable GLVis visualization.");
119  args.Parse();
120  if (!args.Good())
121  {
122  if (myid == 0)
123  {
124  args.PrintUsage(cout);
125  }
126  MPI_Finalize();
127  return 1;
128  }
129  if (myid == 0)
130  {
131  args.PrintOptions(cout);
132  }
133 
134  MFEM_VERIFY(prob >= 0 && prob <=2,
135  "Unrecognized problem type: " << prob);
136 
137  if ( a_coef != 0.0 )
138  {
139  mu_ = 1.0 / a_coef;
140  }
141  if ( freq > 0.0 )
142  {
143  omega_ = 2.0 * M_PI * freq;
144  }
145 
146  exact_sol = check_for_inline_mesh(mesh_file);
147  if (myid == 0 && exact_sol)
148  {
149  cout << "Identified a mesh with known exact solution" << endl;
150  }
151 
153  herm_conv ? ComplexOperator::HERMITIAN : ComplexOperator::BLOCK_SYMMETRIC;
154 
155  // 3. Read the (serial) mesh from the given mesh file on all processors. We
156  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
157  // and volume meshes with the same code.
158  Mesh *mesh = new Mesh(mesh_file, 1, 1);
159  int dim = mesh->Dimension();
160 
161  // 4. Refine the serial mesh on all processors to increase the resolution.
162  for (int l = 0; l < ser_ref_levels; l++)
163  {
164  mesh->UniformRefinement();
165  }
166 
167  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
168  // this mesh further in parallel to increase the resolution. Once the
169  // parallel mesh is defined, the serial mesh can be deleted.
170  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
171  delete mesh;
172  for (int l = 0; l < par_ref_levels; l++)
173  {
174  pmesh->UniformRefinement();
175  }
176 
177  // 6. Define a parallel finite element space on the parallel mesh. Here we
178  // use continuous Lagrange, Nedelec, or Raviart-Thomas finite elements of
179  // the specified order.
180  if (dim == 1 && prob != 0 )
181  {
182  if (myid == 0)
183  {
184  cout << "Switching to problem type 0, H1 basis functions, "
185  << "for 1 dimensional mesh." << endl;
186  }
187  prob = 0;
188  }
189 
190  FiniteElementCollection *fec = NULL;
191  switch (prob)
192  {
193  case 0: fec = new H1_FECollection(order, dim); break;
194  case 1: fec = new ND_FECollection(order, dim); break;
195  case 2: fec = new RT_FECollection(order - 1, dim); break;
196  default: break; // This should be unreachable
197  }
198  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
199  HYPRE_Int size = fespace->GlobalTrueVSize();
200  if (myid == 0)
201  {
202  cout << "Number of finite element unknowns: " << size << endl;
203  }
204 
205  // 7. Determine the list of true (i.e. parallel conforming) essential
206  // boundary dofs. In this example, the boundary conditions are defined
207  // based on the type of mesh and the problem type.
208  Array<int> ess_tdof_list;
209  Array<int> ess_bdr;
210  if (pmesh->bdr_attributes.Size())
211  {
212  ess_bdr.SetSize(pmesh->bdr_attributes.Max());
213  ess_bdr = 1;
214  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
215  }
216 
217  // 8. Set up the parallel linear form b(.) which corresponds to the
218  // right-hand side of the FEM linear system.
219  ParComplexLinearForm b(fespace, conv);
220  b.Vector::operator=(0.0);
221 
222  // 9. Define the solution vector u as a parallel complex finite element grid
223  // function corresponding to fespace. Initialize u with initial guess of
224  // 1+0i or the exact solution if it is known.
225  ParComplexGridFunction u(fespace);
226  ParComplexGridFunction * u_exact = NULL;
227  if (exact_sol) { u_exact = new ParComplexGridFunction(fespace); }
228 
235 
236  ConstantCoefficient zeroCoef(0.0);
237  ConstantCoefficient oneCoef(1.0);
238 
239  Vector zeroVec(dim); zeroVec = 0.0;
240  Vector oneVec(dim); oneVec = 0.0; oneVec[(prob==2)?(dim-1):0] = 1.0;
241  VectorConstantCoefficient zeroVecCoef(zeroVec);
242  VectorConstantCoefficient oneVecCoef(oneVec);
243 
244  u = 0.0;
245  switch (prob)
246  {
247  case 0:
248  if (exact_sol)
249  {
250  u.ProjectBdrCoefficient(u0_r, u0_i, ess_bdr);
251  u_exact->ProjectCoefficient(u0_r, u0_i);
252  }
253  else
254  {
255  u.ProjectBdrCoefficient(oneCoef, zeroCoef, ess_bdr);
256  }
257  break;
258  case 1:
259  if (exact_sol)
260  {
261  u.ProjectBdrCoefficientTangent(u1_r, u1_i, ess_bdr);
262  u_exact->ProjectCoefficient(u1_r, u1_i);
263  }
264  else
265  {
266  u.ProjectBdrCoefficientTangent(oneVecCoef, zeroVecCoef, ess_bdr);
267  }
268  break;
269  case 2:
270  if (exact_sol)
271  {
272  u.ProjectBdrCoefficientNormal(u2_r, u2_i, ess_bdr);
273  u_exact->ProjectCoefficient(u2_r, u2_i);
274  }
275  else
276  {
277  u.ProjectBdrCoefficientNormal(oneVecCoef, zeroVecCoef, ess_bdr);
278  }
279  break;
280  default: break; // This should be unreachable
281  }
282 
283  if (visualization && exact_sol)
284  {
285  char vishost[] = "localhost";
286  int visport = 19916;
287  socketstream sol_sock_r(vishost, visport);
288  socketstream sol_sock_i(vishost, visport);
289  sol_sock_r << "parallel " << num_procs << " " << myid << "\n";
290  sol_sock_i << "parallel " << num_procs << " " << myid << "\n";
291  sol_sock_r.precision(8);
292  sol_sock_i.precision(8);
293  sol_sock_r << "solution\n" << *pmesh << u_exact->real()
294  << "window_title 'Exact: Real Part'" << flush;
295  sol_sock_i << "solution\n" << *pmesh << u_exact->imag()
296  << "window_title 'Exact: Imaginary Part'" << flush;
297  }
298 
299  // 10. Set up the parallel sesquilinear form a(.,.) on the finite element
300  // space corresponding to the damped harmonic oscillator operator of the
301  // appropriate type:
302  //
303  // 0) A scalar H1 field
304  // -Div(a Grad) - omega^2 b + i omega c
305  //
306  // 1) A vector H(Curl) field
307  // Curl(a Curl) - omega^2 b + i omega c
308  //
309  // 2) A vector H(Div) field
310  // -Grad(a Div) - omega^2 b + i omega c
311  //
312  ConstantCoefficient stiffnessCoef(1.0/mu_);
313  ConstantCoefficient massCoef(-omega_ * omega_ * epsilon_);
314  ConstantCoefficient lossCoef(omega_ * sigma_);
315  ConstantCoefficient negMassCoef(omega_ * omega_ * epsilon_);
316 
317  ParSesquilinearForm *a = new ParSesquilinearForm(fespace, conv);
318  switch (prob)
319  {
320  case 0:
321  a->AddDomainIntegrator(new DiffusionIntegrator(stiffnessCoef),
322  NULL);
323  a->AddDomainIntegrator(new MassIntegrator(massCoef),
324  new MassIntegrator(lossCoef));
325  break;
326  case 1:
327  a->AddDomainIntegrator(new CurlCurlIntegrator(stiffnessCoef),
328  NULL);
330  new VectorFEMassIntegrator(lossCoef));
331  break;
332  case 2:
333  a->AddDomainIntegrator(new DivDivIntegrator(stiffnessCoef),
334  NULL);
336  new VectorFEMassIntegrator(lossCoef));
337  break;
338  default: break; // This should be unreachable
339  }
340 
341  // 10a. Set up the parallel bilinear form for the preconditioner
342  // corresponding to the appropriate operator
343  //
344  // 0) A scalar H1 field
345  // -Div(a Grad) - omega^2 b + omega c
346  //
347  // 1) A vector H(Curl) field
348  // Curl(a Curl) + omega^2 b + omega c
349  //
350  // 2) A vector H(Div) field
351  // -Grad(a Div) - omega^2 b + omega c
352  //
353  ParBilinearForm *pcOp = new ParBilinearForm(fespace);
354  switch (prob)
355  {
356  case 0:
357  pcOp->AddDomainIntegrator(new DiffusionIntegrator(stiffnessCoef));
358  pcOp->AddDomainIntegrator(new MassIntegrator(massCoef));
359  pcOp->AddDomainIntegrator(new MassIntegrator(lossCoef));
360  break;
361  case 1:
362  pcOp->AddDomainIntegrator(new CurlCurlIntegrator(stiffnessCoef));
363  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(negMassCoef));
364  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(lossCoef));
365  break;
366  case 2:
367  pcOp->AddDomainIntegrator(new DivDivIntegrator(stiffnessCoef));
368  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(massCoef));
369  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(lossCoef));
370  break;
371  default: break; // This should be unreachable
372  }
373 
374  // 11. Assemble the parallel bilinear form and the corresponding linear
375  // system, applying any necessary transformations such as: parallel
376  // assembly, eliminating boundary conditions, applying conforming
377  // constraints for non-conforming AMR, etc.
378  a->Assemble();
379  pcOp->Assemble();
380 
381  OperatorHandle A;
382  Vector B, U;
383 
384  a->FormLinearSystem(ess_tdof_list, u, b, A, U, B);
385  u = 0.0;
386  U = 0.0;
387 
388  OperatorHandle PCOp;
389  pcOp->FormSystemMatrix(ess_tdof_list, PCOp);
390 
391  if (myid == 0)
392  {
393  ComplexHypreParMatrix * Ahyp =
394  dynamic_cast<ComplexHypreParMatrix*>(A.Ptr());
395 
396  cout << "Size of linear system: "
397  << 2 * Ahyp->real().GetGlobalNumRows() << endl << endl;
398  }
399 
400  // 12. Define and apply a parallel FGMRES solver for AU=B with a block
401  // diagonal preconditioner based on the appropriate multigrid
402  // preconditioner from hypre.
403  {
404  Array<int> blockTrueOffsets;
405  blockTrueOffsets.SetSize(3);
406  blockTrueOffsets[0] = 0;
407  blockTrueOffsets[1] = PCOp.Ptr()->Height();
408  blockTrueOffsets[2] = PCOp.Ptr()->Height();
409  blockTrueOffsets.PartialSum();
410 
411  BlockDiagonalPreconditioner BDP(blockTrueOffsets);
412 
413  Operator * pc_r = NULL;
414  Operator * pc_i = NULL;
415 
416  switch (prob)
417  {
418  case 0:
419  pc_r = new HypreBoomerAMG(*PCOp.As<HypreParMatrix>());
420  break;
421  case 1:
422  pc_r = new HypreAMS(*PCOp.As<HypreParMatrix>(), fespace);
423  break;
424  case 2:
425  if (dim == 2 )
426  {
427  pc_r = new HypreAMS(*PCOp.As<HypreParMatrix>(), fespace);
428  }
429  else
430  {
431  pc_r = new HypreADS(*PCOp.As<HypreParMatrix>(), fespace);
432  }
433  break;
434  default: break; // This should be unreachable
435  }
436  pc_i = new ScaledOperator(pc_r,
437  (conv == ComplexOperator::HERMITIAN) ?
438  -1.0:1.0);
439 
440  BDP.SetDiagonalBlock(0, pc_r);
441  BDP.SetDiagonalBlock(1, pc_i);
442  BDP.owns_blocks = 1;
443 
444  FGMRESSolver fgmres(MPI_COMM_WORLD);
445  fgmres.SetPreconditioner(BDP);
446  fgmres.SetOperator(*A.Ptr());
447  fgmres.SetRelTol(1e-12);
448  fgmres.SetMaxIter(1000);
449  fgmres.SetPrintLevel(1);
450  fgmres.Mult(B, U);
451  }
452  // 13. Recover the parallel grid function corresponding to U. This is the
453  // local finite element solution on each processor.
454  a->RecoverFEMSolution(U, b, u);
455 
456  if (exact_sol)
457  {
458  double err_r = -1.0;
459  double err_i = -1.0;
460 
461  switch (prob)
462  {
463  case 0:
464  err_r = u.real().ComputeL2Error(u0_r);
465  err_i = u.imag().ComputeL2Error(u0_i);
466  break;
467  case 1:
468  err_r = u.real().ComputeL2Error(u1_r);
469  err_i = u.imag().ComputeL2Error(u1_i);
470  break;
471  case 2:
472  err_r = u.real().ComputeL2Error(u2_r);
473  err_i = u.imag().ComputeL2Error(u2_i);
474  break;
475  default: break; // This should be unreachable
476  }
477 
478  if ( myid == 0 )
479  {
480  cout << endl;
481  cout << "|| Re (u_h - u) ||_{L^2} = " << err_r << endl;
482  cout << "|| Im (u_h - u) ||_{L^2} = " << err_i << endl;
483  cout << endl;
484  }
485  }
486 
487  // 14. Save the refined mesh and the solution in parallel. This output can be
488  // viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
489  {
490  ostringstream mesh_name, sol_r_name, sol_i_name;
491  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
492  sol_r_name << "sol_r." << setfill('0') << setw(6) << myid;
493  sol_i_name << "sol_i." << setfill('0') << setw(6) << myid;
494 
495  ofstream mesh_ofs(mesh_name.str().c_str());
496  mesh_ofs.precision(8);
497  pmesh->Print(mesh_ofs);
498 
499  ofstream sol_r_ofs(sol_r_name.str().c_str());
500  ofstream sol_i_ofs(sol_i_name.str().c_str());
501  sol_r_ofs.precision(8);
502  sol_i_ofs.precision(8);
503  u.real().Save(sol_r_ofs);
504  u.imag().Save(sol_i_ofs);
505  }
506 
507  // 15. Send the solution by socket to a GLVis server.
508  if (visualization)
509  {
510  char vishost[] = "localhost";
511  int visport = 19916;
512  socketstream sol_sock_r(vishost, visport);
513  socketstream sol_sock_i(vishost, visport);
514  sol_sock_r << "parallel " << num_procs << " " << myid << "\n";
515  sol_sock_i << "parallel " << num_procs << " " << myid << "\n";
516  sol_sock_r.precision(8);
517  sol_sock_i.precision(8);
518  sol_sock_r << "solution\n" << *pmesh << u.real()
519  << "window_title 'Solution: Real Part'" << flush;
520  sol_sock_i << "solution\n" << *pmesh << u.imag()
521  << "window_title 'Solution: Imaginary Part'" << flush;
522  }
523  if (visualization && exact_sol)
524  {
525  *u_exact -= u;
526 
527  char vishost[] = "localhost";
528  int visport = 19916;
529  socketstream sol_sock_r(vishost, visport);
530  socketstream sol_sock_i(vishost, visport);
531  sol_sock_r << "parallel " << num_procs << " " << myid << "\n";
532  sol_sock_i << "parallel " << num_procs << " " << myid << "\n";
533  sol_sock_r.precision(8);
534  sol_sock_i.precision(8);
535  sol_sock_r << "solution\n" << *pmesh << u_exact->real()
536  << "window_title 'Error: Real Part'" << flush;
537  sol_sock_i << "solution\n" << *pmesh << u_exact->imag()
538  << "window_title 'Error: Imaginary Part'" << flush;
539  }
540  if (visualization)
541  {
542  ParGridFunction u_t(fespace);
543  u_t = u.real();
544  char vishost[] = "localhost";
545  int visport = 19916;
546  socketstream sol_sock(vishost, visport);
547  sol_sock << "parallel " << num_procs << " " << myid << "\n";
548  sol_sock.precision(8);
549  sol_sock << "solution\n" << *pmesh << u_t
550  << "window_title 'Harmonic Solution (t = 0.0 T)'"
551  << "pause\n" << flush;
552  if (myid == 0)
553  cout << "GLVis visualization paused."
554  << " Press space (in the GLVis window) to resume it.\n";
555  int num_frames = 32;
556  int i = 0;
557  while (sol_sock)
558  {
559  double t = (double)(i % num_frames) / num_frames;
560  ostringstream oss;
561  oss << "Harmonic Solution (t = " << t << " T)";
562 
563  add(cos( 2.0 * M_PI * t), u.real(),
564  sin(-2.0 * M_PI * t), u.imag(), u_t);
565  sol_sock << "parallel " << num_procs << " " << myid << "\n";
566  sol_sock << "solution\n" << *pmesh << u_t
567  << "window_title '" << oss.str() << "'" << flush;
568  i++;
569  }
570  }
571 
572  // 16. Free the used memory.
573  delete a;
574  delete u_exact;
575  delete pcOp;
576  delete fespace;
577  delete fec;
578  delete pmesh;
579 
580  MPI_Finalize();
581 
582  return 0;
583 }
584 
585 bool check_for_inline_mesh(const char * mesh_file)
586 {
587  string file(mesh_file);
588  size_t p0 = file.find_last_of("/");
589  string s0 = file.substr((p0==string::npos)?0:(p0+1),7);
590  return s0 == "inline-";
591 }
592 
593 complex<double> u0_exact(const Vector &x)
594 {
595  int dim = x.Size();
596  complex<double> i(0.0, 1.0);
597  complex<double> alpha = (epsilon_ * omega_ - i * sigma_);
598  complex<double> kappa = std::sqrt(mu_ * omega_* alpha);
599  return std::exp(-i * kappa * x[dim - 1]);
600 }
601 
602 double u0_real_exact(const Vector &x)
603 {
604  return u0_exact(x).real();
605 }
606 
607 double u0_imag_exact(const Vector &x)
608 {
609  return u0_exact(x).imag();
610 }
611 
612 void u1_real_exact(const Vector &x, Vector &v)
613 {
614  int dim = x.Size();
615  v.SetSize(dim); v = 0.0; v[0] = u0_real_exact(x);
616 }
617 
618 void u1_imag_exact(const Vector &x, Vector &v)
619 {
620  int dim = x.Size();
621  v.SetSize(dim); v = 0.0; v[0] = u0_imag_exact(x);
622 }
623 
624 void u2_real_exact(const Vector &x, Vector &v)
625 {
626  int dim = x.Size();
627  v.SetSize(dim); v = 0.0; v[dim-1] = u0_real_exact(x);
628 }
629 
630 void u2_imag_exact(const Vector &x, Vector &v)
631 {
632  int dim = x.Size();
633  v.SetSize(dim); v = 0.0; v[dim-1] = u0_imag_exact(x);
634 }
int Size() const
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:770
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:96
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1014
The Auxiliary-space Divergence Solver in hypre.
Definition: hypre.hpp:1053
Subclass constant coefficient.
Definition: coefficient.hpp:67
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
void u2_imag_exact(const Vector &, Vector &)
Definition: ex22.cpp:575
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
Integrator for (curl u, curl v) for Nedelec elements.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
FGMRES method.
Definition: solvers.hpp:203
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:774
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
ParGridFunction & imag()
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
Abstract parallel finite element space.
Definition: pfespace.hpp:28
complex< double > u0_exact(const Vector &x)
Definition: ex22.cpp:538
int main(int argc, char *argv[])
Definition: ex1.cpp:62
(Q div u, div v) for RT elements
void u1_imag_exact(const Vector &, Vector &)
Definition: ex22.cpp:563
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:238
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:421
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
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
double u0_imag_exact(const Vector &)
Definition: ex22.cpp:552
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
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.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4102
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:191
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
void SetRelTol(double rtol)
Definition: solvers.hpp:63
Scaled Operator B: x -&gt; a A(x).
Definition: operator.hpp:544
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)
Definition: optparser.hpp:76
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
void PartialSum()
Partial Sum.
Definition: array.cpp:103
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
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.
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL) const
Definition: pgridfunc.hpp:250
int dim
Definition: ex24.cpp:43
Specialization of the ComplexOperator built from a pair of HypreParMatrices.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
const double alpha
Definition: ex15.cpp:336
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:251
double kappa
Definition: ex3.cpp:54
Vector data type.
Definition: vector.hpp:48
void u1_real_exact(const Vector &, Vector &)
Definition: ex22.cpp:557
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
ParGridFunction & real()
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Abstract operator.
Definition: operator.hpp:24
bool check_for_inline_mesh(const char *mesh_file)
Definition: ex22.cpp:530
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
double freq
Definition: ex3.cpp:54
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual HypreParMatrix & real()
Real or imaginary part accessor methods.
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
void u2_real_exact(const Vector &, Vector &)
Definition: ex22.cpp:569
bool Good() const
Definition: optparser.hpp:125
double u0_real_exact(const Vector &)
Definition: ex22.cpp:547