MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex31p.cpp
Go to the documentation of this file.
1 // MFEM Example 31 - Parallel Version
2 //
3 // Compile with: make ex31p
4 //
5 // Sample runs: mpirun -np 4 ex31p -m ../data/hexagon.mesh -o 2
6 // mpirun -np 4 ex31p -m ../data/star.mesh
7 // mpirun -np 4 ex31p -m ../data/square-disc.mesh -o 2
8 // mpirun -np 4 ex31p -m ../data/fichera.mesh -o 3 -rs 1 -rp 0
9 // mpirun -np 4 ex31p -m ../data/square-disc-nurbs.mesh -o 3
10 // mpirun -np 4 ex31p -m ../data/amr-quad.mesh -o 2 -rs 1
11 // mpirun -np 4 ex31p -m ../data/amr-hex.mesh -rs 1
12 //
13 // Description: This example code solves a simple electromagnetic diffusion
14 // problem corresponding to the second order definite Maxwell
15 // equation curl curl E + sigma E = f with boundary condition
16 // E x n = <given tangential field>. In this example sigma is an
17 // anisotropic 3x3 tensor. Here, we use a given exact solution E
18 // and compute the corresponding r.h.s. f. We discretize with
19 // Nedelec finite elements in 1D, 2D, or 3D.
20 //
21 // The example demonstrates the use of restricted H(curl) finite
22 // element spaces with the curl-curl and the (vector finite
23 // element) mass bilinear form, as well as the computation of
24 // discretization error when the exact solution is known. These
25 // restricted spaces allow the solution of 1D or 2D
26 // electromagnetic problems which involve 3D field vectors. Such
27 // problems arise in plasma physics and crystallography.
28 //
29 // We recommend viewing example 3 before viewing this example.
30 
31 #include "mfem.hpp"
32 #include <fstream>
33 #include <iostream>
34 
35 using namespace std;
36 using namespace mfem;
37 
38 // Exact solution, E, and r.h.s., f. See below for implementation.
39 void E_exact(const Vector &, Vector &);
40 void CurlE_exact(const Vector &, Vector &);
41 void f_exact(const Vector &, Vector &);
42 double freq = 1.0, kappa;
43 int dim;
44 
45 int main(int argc, char *argv[])
46 {
47  // 1. Initialize MPI.
48  MPI_Session mpi;
49  int num_procs = mpi.WorldSize();
50  int myid = mpi.WorldRank();
51 
52  // 2. Parse command-line options.
53  const char *mesh_file = "../data/inline-quad.mesh";
54  int ser_ref_levels = 2;
55  int par_ref_levels = 1;
56  int order = 1;
57  bool use_ams = true;
58  bool visualization = 1;
59 
60  OptionsParser args(argc, argv);
61  args.AddOption(&mesh_file, "-m", "--mesh",
62  "Mesh file to use.");
63  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
64  "Number of times to refine the mesh uniformly in serial.");
65  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
66  "Number of times to refine the mesh uniformly in parallel.");
67  args.AddOption(&order, "-o", "--order",
68  "Finite element order (polynomial degree).");
69  args.AddOption(&freq, "-f", "--frequency", "Set the frequency for the exact"
70  " solution.");
71  args.AddOption(&use_ams, "-ams", "--hypre-ams", "-slu",
72  "--superlu", "Use AMS or SuperLU solver.");
73  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
74  "--no-visualization",
75  "Enable or disable GLVis visualization.");
76  args.ParseCheck();
77 
78  kappa = freq * M_PI;
79 
80  // 3. Read the (serial) mesh from the given mesh file on all processors. We
81  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
82  // and volume meshes with the same code.
83  Mesh *mesh = new Mesh(mesh_file, 1, 1);
84  dim = mesh->Dimension();
85 
86  // 4. Refine the serial mesh on all processors to increase the resolution. In
87  // this example we do 'ref_levels' of uniform refinement (2 by default, or
88  // specified on the command line with -rs).
89  for (int lev = 0; lev < ser_ref_levels; lev++)
90  {
91  mesh->UniformRefinement();
92  }
93 
94  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
95  // this mesh further in parallel to increase the resolution (1 time by
96  // default, or specified on the command line with -rp). Once the parallel
97  // mesh is defined, the serial mesh can be deleted.
98  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
99  delete mesh;
100  for (int lev = 0; lev < par_ref_levels; lev++)
101  {
102  pmesh.UniformRefinement();
103  }
104 
105  // 6. Define a parallel finite element space on the parallel mesh. Here we
106  // use the Nedelec finite elements of the specified order.
107  FiniteElementCollection *fec = NULL;
108  if (dim == 1)
109  {
110  fec = new ND_R1D_FECollection(order, dim);
111  }
112  else if (dim == 2)
113  {
114  fec = new ND_R2D_FECollection(order, dim);
115  }
116  else
117  {
118  fec = new ND_FECollection(order, dim);
119  }
120  ParFiniteElementSpace fespace(&pmesh, fec);
121  HYPRE_Int size = fespace.GlobalTrueVSize();
122  if (mpi.Root()) { cout << "Number of H(Curl) unknowns: " << size << endl; }
123 
124  // 7. Determine the list of true (i.e. parallel conforming) essential
125  // boundary dofs. In this example, the boundary conditions are defined
126  // by marking all the boundary attributes from the mesh as essential
127  // (Dirichlet) and converting them to a list of true dofs.
128  Array<int> ess_tdof_list;
129  if (pmesh.bdr_attributes.Size())
130  {
131  Array<int> ess_bdr(pmesh.bdr_attributes.Max());
132  ess_bdr = 1;
133  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
134  }
135 
136  // 8. Set up the parallel linear form b(.) which corresponds to the
137  // right-hand side of the FEM linear system, which in this case is
138  // (f,phi_i) where f is given by the function f_exact and phi_i are the
139  // basis functions in the finite element fespace.
141  ParLinearForm b(&fespace);
143  b.Assemble();
144 
145  // 9. Define the solution vector x as a parallel finite element grid function
146  // corresponding to fespace. Initialize x by projecting the exact
147  // solution. Note that only values from the boundary edges will be used
148  // when eliminating the non-homogeneous boundary condition to modify the
149  // r.h.s. vector b.
150  ParGridFunction sol(&fespace);
153  sol.ProjectCoefficient(E);
154 
155  // 10. Set up the parallel bilinear form corresponding to the EM diffusion
156  // operator curl muinv curl + sigma I, by adding the curl-curl and the
157  // mass domain integrators.
158  DenseMatrix sigmaMat(3);
159  sigmaMat(0,0) = 2.0; sigmaMat(1,1) = 2.0; sigmaMat(2,2) = 2.0;
160  sigmaMat(0,2) = 0.0; sigmaMat(2,0) = 0.0;
161  sigmaMat(0,1) = M_SQRT1_2; sigmaMat(1,0) = M_SQRT1_2; // 1/sqrt(2) in cmath
162  sigmaMat(1,2) = M_SQRT1_2; sigmaMat(2,1) = M_SQRT1_2;
163 
164  ConstantCoefficient muinv(1.0);
166  ParBilinearForm a(&fespace);
169 
170  // 11. Assemble the parallel bilinear form and the corresponding linear
171  // system, applying any necessary transformations such as: parallel
172  // assembly, eliminating boundary conditions, applying conforming
173  // constraints for non-conforming AMR, etc.
174  a.Assemble();
175 
176  OperatorPtr A;
177  Vector B, X;
178 
179  a.FormLinearSystem(ess_tdof_list, sol, b, A, X, B);
180 
181  // 12. Solve the system AX=B using PCG with the AMS preconditioner from hypre
182  if (use_ams)
183  {
184  if (mpi.Root())
185  {
186  cout << "Size of linear system: "
187  << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
188  }
189 
190  HypreAMS ams(*A.As<HypreParMatrix>(), &fespace);
191 
192  HyprePCG pcg(*A.As<HypreParMatrix>());
193  pcg.SetTol(1e-12);
194  pcg.SetMaxIter(1000);
195  pcg.SetPrintLevel(2);
196  pcg.SetPreconditioner(ams);
197  pcg.Mult(B, X);
198  }
199  else
200 #ifdef MFEM_USE_SUPERLU
201  {
202  if (mpi.Root())
203  {
204  cout << "Size of linear system: "
205  << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
206  }
207 
208  SuperLURowLocMatrix A_SuperLU(*A.As<HypreParMatrix>());
209  SuperLUSolver AInv(MPI_COMM_WORLD);
210  AInv.SetOperator(A_SuperLU);
211  AInv.Mult(B,X);
212  }
213 #else
214  {
215  if (mpi.Root()) { cout << "No solvers available." << endl; }
216  return 1;
217  }
218 #endif
219 
220  // 13. Recover the parallel grid function corresponding to X. This is the
221  // local finite element solution on each processor.
222  a.RecoverFEMSolution(X, b, sol);
223 
224  // 14. Compute and print the H(Curl) norm of the error.
225  {
226  double error = sol.ComputeHCurlError(&E, &CurlE);
227  if (mpi.Root())
228  {
229  cout << "\n|| E_h - E ||_{H(Curl)} = " << error << '\n' << endl;
230  }
231  }
232 
233 
234  // 15. Save the refined mesh and the solution in parallel. This output can
235  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
236  {
237  ostringstream mesh_name, sol_name;
238  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
239  sol_name << "sol." << setfill('0') << setw(6) << myid;
240 
241  ofstream mesh_ofs(mesh_name.str().c_str());
242  mesh_ofs.precision(8);
243  pmesh.Print(mesh_ofs);
244 
245  ofstream sol_ofs(sol_name.str().c_str());
246  sol_ofs.precision(8);
247  sol.Save(sol_ofs);
248  }
249 
250  // 16. Send the solution by socket to a GLVis server.
251  if (visualization)
252  {
253  char vishost[] = "localhost";
254  int visport = 19916;
255 
256  VectorGridFunctionCoefficient solCoef(&sol);
257  CurlGridFunctionCoefficient dsolCoef(&sol);
258 
259  if (dim ==1)
260  {
261  socketstream x_sock(vishost, visport);
262  socketstream y_sock(vishost, visport);
263  socketstream z_sock(vishost, visport);
264  socketstream dy_sock(vishost, visport);
265  socketstream dz_sock(vishost, visport);
266  x_sock.precision(8);
267  y_sock.precision(8);
268  z_sock.precision(8);
269  dy_sock.precision(8);
270  dz_sock.precision(8);
271 
272  Vector xVec(3); xVec = 0.0; xVec(0) = 1;
273  Vector yVec(3); yVec = 0.0; yVec(1) = 1;
274  Vector zVec(3); zVec = 0.0; zVec(2) = 1;
275  VectorConstantCoefficient xVecCoef(xVec);
276  VectorConstantCoefficient yVecCoef(yVec);
277  VectorConstantCoefficient zVecCoef(zVec);
278 
279  H1_FECollection fec_h1(order, dim);
280  L2_FECollection fec_l2(order-1, dim);
281 
282  ParFiniteElementSpace fes_h1(&pmesh, &fec_h1);
283  ParFiniteElementSpace fes_l2(&pmesh, &fec_l2);
284 
285  ParGridFunction xComp(&fes_l2);
286  ParGridFunction yComp(&fes_h1);
287  ParGridFunction zComp(&fes_h1);
288 
289  ParGridFunction dyComp(&fes_l2);
290  ParGridFunction dzComp(&fes_l2);
291 
292  InnerProductCoefficient xCoef(xVecCoef, solCoef);
293  InnerProductCoefficient yCoef(yVecCoef, solCoef);
294  InnerProductCoefficient zCoef(zVecCoef, solCoef);
295 
296  xComp.ProjectCoefficient(xCoef);
297  yComp.ProjectCoefficient(yCoef);
298  zComp.ProjectCoefficient(zCoef);
299 
300  x_sock << "parallel " << num_procs << " " << myid << "\n"
301  << "solution\n" << pmesh << xComp << flush
302  << "window_title 'X component'" << endl;
303  y_sock << "parallel " << num_procs << " " << myid << "\n"
304  << "solution\n" << pmesh << yComp << flush
305  << "window_geometry 403 0 400 350 "
306  << "window_title 'Y component'" << endl;
307  z_sock << "parallel " << num_procs << " " << myid << "\n"
308  << "solution\n" << pmesh << zComp << flush
309  << "window_geometry 806 0 400 350 "
310  << "window_title 'Z component'" << endl;
311 
312  InnerProductCoefficient dyCoef(yVecCoef, dsolCoef);
313  InnerProductCoefficient dzCoef(zVecCoef, dsolCoef);
314 
315  dyComp.ProjectCoefficient(dyCoef);
316  dzComp.ProjectCoefficient(dzCoef);
317 
318  dy_sock << "parallel " << num_procs << " " << myid << "\n"
319  << "solution\n" << pmesh << dyComp << flush
320  << "window_geometry 403 375 400 350 "
321  << "window_title 'Y component of Curl'" << endl;
322  dz_sock << "parallel " << num_procs << " " << myid << "\n"
323  << "solution\n" << pmesh << dzComp << flush
324  << "window_geometry 806 375 400 350 "
325  << "window_title 'Z component of Curl'" << endl;
326  }
327  else if (dim == 2)
328  {
329  socketstream xy_sock(vishost, visport);
330  socketstream z_sock(vishost, visport);
331  socketstream dxy_sock(vishost, visport);
332  socketstream dz_sock(vishost, visport);
333 
334  DenseMatrix xyMat(2,3); xyMat = 0.0;
335  xyMat(0,0) = 1.0; xyMat(1,1) = 1.0;
336  MatrixConstantCoefficient xyMatCoef(xyMat);
337  Vector zVec(3); zVec = 0.0; zVec(2) = 1;
338  VectorConstantCoefficient zVecCoef(zVec);
339 
340  MatrixVectorProductCoefficient xyCoef(xyMatCoef, solCoef);
341  InnerProductCoefficient zCoef(zVecCoef, solCoef);
342 
343  H1_FECollection fec_h1(order, dim);
344  ND_FECollection fec_nd(order, dim);
345  RT_FECollection fec_rt(order-1, dim);
346  L2_FECollection fec_l2(order-1, dim);
347 
348  ParFiniteElementSpace fes_h1(&pmesh, &fec_h1);
349  ParFiniteElementSpace fes_nd(&pmesh, &fec_nd);
350  ParFiniteElementSpace fes_rt(&pmesh, &fec_rt);
351  ParFiniteElementSpace fes_l2(&pmesh, &fec_l2);
352 
353  ParGridFunction xyComp(&fes_nd);
354  ParGridFunction zComp(&fes_h1);
355 
356  ParGridFunction dxyComp(&fes_rt);
357  ParGridFunction dzComp(&fes_l2);
358 
359  xyComp.ProjectCoefficient(xyCoef);
360  zComp.ProjectCoefficient(zCoef);
361 
362  xy_sock << "parallel " << num_procs << " " << myid << "\n";
363  xy_sock.precision(8);
364  xy_sock << "solution\n" << pmesh << xyComp
365  << "window_title 'XY components'\n" << flush;
366  z_sock << "parallel " << num_procs << " " << myid << "\n"
367  << "solution\n" << pmesh << zComp << flush
368  << "window_geometry 403 0 400 350 "
369  << "window_title 'Z component'" << endl;
370 
371  MatrixVectorProductCoefficient dxyCoef(xyMatCoef, dsolCoef);
372  InnerProductCoefficient dzCoef(zVecCoef, dsolCoef);
373 
374  dxyComp.ProjectCoefficient(dxyCoef);
375  dzComp.ProjectCoefficient(dzCoef);
376 
377  dxy_sock << "parallel " << num_procs << " " << myid << "\n"
378  << "solution\n" << pmesh << dxyComp << flush
379  << "window_geometry 0 375 400 350 "
380  << "window_title 'XY components of Curl'" << endl;
381  dz_sock << "parallel " << num_procs << " " << myid << "\n"
382  << "solution\n" << pmesh << dzComp << flush
383  << "window_geometry 403 375 400 350 "
384  << "window_title 'Z component of Curl'" << endl;
385  }
386  else
387  {
388  socketstream sol_sock(vishost, visport);
389  socketstream dsol_sock(vishost, visport);
390 
391  RT_FECollection fec_rt(order-1, dim);
392 
393  ParFiniteElementSpace fes_rt(&pmesh, &fec_rt);
394 
395  ParGridFunction dsol(&fes_rt);
396 
397  dsol.ProjectCoefficient(dsolCoef);
398 
399  sol_sock << "parallel " << num_procs << " " << myid << "\n";
400  sol_sock.precision(8);
401  sol_sock << "solution\n" << pmesh << sol
402  << "window_title 'Solution'" << flush << endl;
403  dsol_sock << "parallel " << num_procs << " " << myid << "\n"
404  << "solution\n" << pmesh << dsol << flush
405  << "window_geometry 0 375 400 350 "
406  << "window_title 'Curl of solution'" << endl;
407  }
408  }
409 
410  // 17. Free the used memory.
411  delete fec;
412 
413  return 0;
414 }
415 
416 
417 void E_exact(const Vector &x, Vector &E)
418 {
419  if (dim == 1)
420  {
421  E(0) = 1.1 * sin(kappa * x(0) + 0.0 * M_PI);
422  E(1) = 1.2 * sin(kappa * x(0) + 0.4 * M_PI);
423  E(2) = 1.3 * sin(kappa * x(0) + 0.9 * M_PI);
424  }
425  else if (dim == 2)
426  {
427  E(0) = 1.1 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
428  E(1) = 1.2 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
429  E(2) = 1.3 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
430  }
431  else
432  {
433  E(0) = 1.1 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
434  E(1) = 1.2 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
435  E(2) = 1.3 * sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
436  E *= cos(kappa * x(2));
437  }
438 }
439 
440 void CurlE_exact(const Vector &x, Vector &dE)
441 {
442  if (dim == 1)
443  {
444  double c4 = cos(kappa * x(0) + 0.4 * M_PI);
445  double c9 = cos(kappa * x(0) + 0.9 * M_PI);
446 
447  dE(0) = 0.0;
448  dE(1) = -1.3 * c9;
449  dE(2) = 1.2 * c4;
450  dE *= kappa;
451  }
452  else if (dim == 2)
453  {
454  double c0 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
455  double c4 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
456  double c9 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
457 
458  dE(0) = 1.3 * c9;
459  dE(1) = -1.3 * c9;
460  dE(2) = 1.2 * c4 - 1.1 * c0;
461  dE *= kappa * M_SQRT1_2;
462  }
463  else
464  {
465  double s0 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
466  double c0 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
467  double s4 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
468  double c4 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
469  double c9 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
470  double sk = sin(kappa * x(2));
471  double ck = cos(kappa * x(2));
472 
473  dE(0) = 1.2 * s4 * sk + 1.3 * M_SQRT1_2 * c9 * ck;
474  dE(1) = -1.1 * s0 * sk - 1.3 * M_SQRT1_2 * c9 * ck;
475  dE(2) = -M_SQRT1_2 * (1.1 * c0 - 1.2 * c4) * ck;
476  dE *= kappa;
477  }
478 }
479 
480 void f_exact(const Vector &x, Vector &f)
481 {
482  if (dim == 1)
483  {
484  double s0 = sin(kappa * x(0) + 0.0 * M_PI);
485  double s4 = sin(kappa * x(0) + 0.4 * M_PI);
486  double s9 = sin(kappa * x(0) + 0.9 * M_PI);
487 
488  f(0) = 2.2 * s0 + 1.2 * M_SQRT1_2 * s4;
489  f(1) = 1.2 * (2.0 + kappa * kappa) * s4 +
490  M_SQRT1_2 * (1.1 * s0 + 1.3 * s9);
491  f(2) = 1.3 * (2.0 + kappa * kappa) * s9 + 1.2 * M_SQRT1_2 * s4;
492  }
493  else if (dim == 2)
494  {
495  double s0 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
496  double s4 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
497  double s9 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
498 
499  f(0) = 0.55 * (4.0 + kappa * kappa) * s0 +
500  0.6 * (M_SQRT2 - kappa * kappa) * s4;
501  f(1) = 0.55 * (M_SQRT2 - kappa * kappa) * s0 +
502  0.6 * (4.0 + kappa * kappa) * s4 +
503  0.65 * M_SQRT2 * s9;
504  f(2) = 0.6 * M_SQRT2 * s4 + 1.3 * (2.0 + kappa * kappa) * s9;
505  }
506  else
507  {
508  double s0 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
509  double c0 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.0 * M_PI);
510  double s4 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
511  double c4 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.4 * M_PI);
512  double s9 = sin(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
513  double c9 = cos(kappa * M_SQRT1_2 * (x(0) + x(1)) + 0.9 * M_PI);
514  double sk = sin(kappa * x(2));
515  double ck = cos(kappa * x(2));
516 
517  f(0) = 0.55 * (4.0 + 3.0 * kappa * kappa) * s0 * ck +
518  0.6 * (M_SQRT2 - kappa * kappa) * s4 * ck -
519  0.65 * M_SQRT2 * kappa * kappa * c9 * sk;
520 
521  f(1) = 0.55 * (M_SQRT2 - kappa * kappa) * s0 * ck +
522  0.6 * (4.0 + 3.0 * kappa * kappa) * s4 * ck +
523  0.65 * M_SQRT2 * s9 * ck -
524  0.65 * M_SQRT2 * kappa * kappa * c9 * sk;
525 
526  f(2) = 0.6 * M_SQRT2 * s4 * ck -
527  M_SQRT2 * kappa * kappa * (0.55 * c0 + 0.6 * c4) * sk
528  + 1.3 * (2.0 + kappa * kappa) * s9 * ck;
529  }
530 }
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
A matrix coefficient that is constant in space and time.
void SetTol(double tol)
Definition: hypre.cpp:3775
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
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Definition: superlu.cpp:474
The Auxiliary-space Maxwell Solver in hypre.
Definition: hypre.hpp:1596
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
Scalar coefficient defined as the inner product of two vector coefficients.
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 1D.
Definition: fe_coll.hpp:460
Integrator for (curl u, curl v) for Nedelec elements.
Vector coefficient that is constant in space and time.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
Arbitrary order 3D H(curl)-conforming Nedelec finite elements in 2D.
Definition: fe_coll.hpp:514
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:525
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:3795
void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: superlu.cpp:580
Class for parallel linear form.
Definition: plinearform.hpp:26
A simple convenience class based on the Mpi singleton class above. Preserved for backward compatibili...
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
constexpr int visport
Vector coefficient defined as the Curl of a vector GridFunction.
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.
bool Root() const
Return true if WorldRank() == 0.
int Dimension() const
Definition: mesh.hpp:999
void SetMaxIter(int max_iter)
Definition: hypre.cpp:3785
FDualNumber< tbase > cos(const FDualNumber< tbase > &f)
cos([dual number])
Definition: fdual.hpp:471
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:337
A general vector function coefficient.
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
int WorldRank() const
Return MPI_COMM_WORLD&#39;s rank.
FDualNumber< tbase > sin(const FDualNumber< tbase > &f)
sin([dual number])
Definition: fdual.hpp:572
PCG solver in hypre.
Definition: hypre.hpp:1116
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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 f_exact(const Vector &, Vector &)
Definition: ex3.cpp:258
Vector coefficient defined as a product of a matrix coefficient and a vector coefficient.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
int dim
Definition: ex24.cpp:53
double freq
Definition: ex24.cpp:54
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:3800
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:242
Class for parallel bilinear form.
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:266
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:411
Vector data type.
Definition: vector.hpp:60
Vector coefficient defined by a vector GridFunction.
void CurlE_exact(const Vector &, Vector &)
Definition: ex31.cpp:375
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4684
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:3823
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
void ParseCheck(std::ostream &out=mfem::out)
Definition: optparser.cpp:252
virtual double ComputeHCurlError(VectorCoefficient *exsol, VectorCoefficient *excurl, const IntegrationRule *irs[]=NULL) const
Returns the error measured H(curl)-norm for ND elements.
Definition: pgridfunc.hpp:362
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:284
double f(const Vector &p)
double sigma(const Vector &x)
Definition: maxwell.cpp:102