MFEM  v4.6.0
Finite element discretization library
ex34.cpp
Go to the documentation of this file.
1 // MFEM Example 34
2 //
3 // Compile with: make ex34
4 //
5 // Sample runs: ex34 -o 2
6 // ex34 -o 2 -pa -hex
7 //
8 // Device sample runs:
9 // ex34 -o 2 -pa -hex -d cuda
10 // ex34 -o 2 -no-pa -d cuda
11 //
12 // Description: This example code solves a simple magnetostatic problem
13 // curl curl A = J where the current density J is computed on a
14 // subset of the domain as J = -sigma grad phi. We discretize the
15 // vector potential with Nedelec finite elements, the scalar
16 // potential with Lagrange finite elements, and the current
17 // density with Raviart-Thomas finite elements.
18 //
19 // The example demonstrates the use of a SubMesh to compute the
20 // scalar potential and its associated current density which is
21 // then transferred to the original mesh and used as a source
22 // function.
23 //
24 // Note that this example takes certain liberties with the
25 // current density which is not necessarily divergence free
26 // as it should be. This was done to focus on the use of the
27 // SubMesh to transfer information between a full mesh and a
28 // sub-domain. A more rigorous implementation might employ an
29 // H(div) saddle point solver to obtain a divergence free J on
30 // the SubMesh. It would then also need to ensure that the r.h.s.
31 // of curl curl A = J does in fact lie in the range of the weak
32 // curl operator by performing a divergence cleaning procedure
33 // before the solve. After divergence cleaning the delta
34 // parameter would probably not be needed.
35 //
36 // This example is designed to make use of a specific mesh which
37 // has a known configuration of elements and boundary attributes.
38 // Other meshes could be used but extra care would be required to
39 // properly define the SubMesh and the necessary boundaries.
40 //
41 // We recommend viewing examples 1 and 3 before viewing this
42 // example.
43 
44 #include "mfem.hpp"
45 #include <fstream>
46 #include <iostream>
47 
48 using namespace std;
49 using namespace mfem;
50 
51 static bool pa_ = false;
52 static bool algebraic_ceed_ = false;
53 
54 void ComputeCurrentDensityOnSubMesh(int order,
55  const Array<int> &phi0_attr,
56  const Array<int> &phi1_attr,
57  const Array<int> &jn_zero_attr,
58  GridFunction &j_cond);
59 
60 int main(int argc, char *argv[])
61 {
62  // 1. Parse command-line options.
63  const char *mesh_file = "../data/fichera-mixed.mesh";
64  Array<int> cond_attr;
65  Array<int> submesh_elems;
66  Array<int> sym_plane_attr;
67  Array<int> phi0_attr;
68  Array<int> phi1_attr;
69  Array<int> jn_zero_attr;
70  int ref_levels = 1;
71  int order = 1;
72  double delta_const = 1e-6;
73  bool mixed = true;
74  bool static_cond = false;
75  const char *device_config = "cpu";
76  bool visualization = true;
77 
78  OptionsParser args(argc, argv);
79  args.AddOption(&ref_levels, "-r", "--refine",
80  "Number of times to refine the mesh uniformly.");
81  args.AddOption(&order, "-o", "--order",
82  "Finite element order (polynomial degree).");
83  args.AddOption(&delta_const, "-mc", "--magnetic-cond",
84  "Magnetic Conductivity");
85  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
86  "--no-static-condensation", "Enable static condensation.");
87  args.AddOption(&mixed, "-mixed", "--mixed-mesh", "-hex",
88  "--hex-mesh", "Mixed mesh of hexahedral mesh.");
89  args.AddOption(&pa_, "-pa", "--partial-assembly", "-no-pa",
90  "--no-partial-assembly", "Enable Partial Assembly.");
91  args.AddOption(&device_config, "-d", "--device",
92  "Device configuration string, see Device::Configure().");
93 #ifdef MFEM_USE_CEED
94  args.AddOption(&algebraic_ceed_, "-a", "--algebraic", "-no-a", "--no-algebraic",
95  "Use algebraic Ceed solver");
96 #endif
97  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
98  "--no-visualization",
99  "Enable or disable GLVis visualization.");
100 
101  args.Parse();
102  if (!args.Good())
103  {
104  args.PrintUsage(cout);
105  return 1;
106  }
107  args.PrintOptions(cout);
108 
109  if (!mixed || pa_)
110  {
111  mesh_file = "../data/fichera.mesh";
112  }
113 
114  if (submesh_elems.Size() == 0)
115  {
116  if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0)
117  {
118  submesh_elems.SetSize(5);
119  submesh_elems[0] = 0;
120  submesh_elems[1] = 2;
121  submesh_elems[2] = 3;
122  submesh_elems[3] = 4;
123  submesh_elems[4] = 9;
124  }
125  else if (strcmp(mesh_file, "../data/fichera.mesh") == 0)
126  {
127  submesh_elems.SetSize(7);
128  submesh_elems[0] = 10;
129  submesh_elems[1] = 14;
130  submesh_elems[2] = 34;
131  submesh_elems[3] = 36;
132  submesh_elems[4] = 37;
133  submesh_elems[5] = 38;
134  submesh_elems[6] = 39;
135  }
136  }
137  if (sym_plane_attr.Size() == 0)
138  {
139  if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
140  strcmp(mesh_file, "../data/fichera.mesh") == 0)
141  {
142  sym_plane_attr.SetSize(8);
143  sym_plane_attr[0] = 9;
144  sym_plane_attr[1] = 10;
145  sym_plane_attr[2] = 11;
146  sym_plane_attr[3] = 12;
147  sym_plane_attr[4] = 13;
148  sym_plane_attr[5] = 14;
149  sym_plane_attr[6] = 15;
150  sym_plane_attr[7] = 16;
151  }
152  }
153  if (phi0_attr.Size() == 0)
154  {
155  if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
156  strcmp(mesh_file, "../data/fichera.mesh") == 0)
157  {
158  phi0_attr.Append(2);
159  }
160  }
161  if (phi1_attr.Size() == 0)
162  {
163  if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
164  strcmp(mesh_file, "../data/fichera.mesh") == 0)
165  {
166  phi1_attr.Append(23);
167  }
168  }
169  if (jn_zero_attr.Size() == 0)
170  {
171  if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
172  strcmp(mesh_file, "../data/fichera.mesh") == 0)
173  {
174  jn_zero_attr.Append(25);
175  }
176  for (int i=0; i<sym_plane_attr.Size(); i++)
177  {
178  jn_zero_attr.Append(sym_plane_attr[i]);
179  }
180  }
181 
182  // 2. Enable hardware devices such as GPUs, and programming models such as
183  // CUDA, OCCA, RAJA and OpenMP based on command line options.
184  Device device(device_config);
185  device.Print();
186 
187  // 3. Read the (serial) mesh from the given mesh file on all processors. We
188  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
189  // and volume meshes with the same code.
190  Mesh mesh(mesh_file, 1, 1);
191  int dim = mesh.Dimension();
192 
193  if (!mixed || pa_)
194  {
195  mesh.UniformRefinement();
196 
197  if (ref_levels > 0)
198  {
199  ref_levels--;
200  }
201  }
202 
203  int submesh_attr = -1;
204  if (cond_attr.Size() == 0 && submesh_elems.Size() > 0)
205  {
206  int max_attr = mesh.attributes.Max();
207  submesh_attr = max_attr + 1;
208 
209  for (int i=0; i<submesh_elems.Size(); i++)
210  {
211  mesh.SetAttribute(submesh_elems[i], submesh_attr);
212  }
213  mesh.SetAttributes();
214 
215  if (cond_attr.Size() == 0)
216  {
217  cond_attr.Append(submesh_attr);
218  }
219  }
220 
221  // 4. Refine the serial mesh on all processors to increase the resolution. In
222  // this example we do 'ref_levels' of uniform refinement.
223  {
224  for (int l = 0; l < ref_levels; l++)
225  {
226  mesh.UniformRefinement();
227  }
228  }
229 
230  // 5b. Extract a submesh covering a portion of the domain
231  SubMesh mesh_cond(SubMesh::CreateFromDomain(mesh, cond_attr));
232 
233  // 6. Define a suitable finite element space on the SubMesh and compute
234  // the current density as an H(div) field.
235  RT_FECollection fec_cond_rt(order - 1, dim);
236  FiniteElementSpace fes_cond_rt(&mesh_cond, &fec_cond_rt);
237  GridFunction j_cond(&fes_cond_rt);
238 
239  ComputeCurrentDensityOnSubMesh(order, phi0_attr, phi1_attr, jn_zero_attr,
240  j_cond);
241 
242  // 6a. Save the SubMesh and associated current density in parallel. This
243  // output can be viewed later using GLVis:
244  // "glvis -np <np> -m cond_mesh -g cond_j"
245  {
246  ostringstream mesh_name, cond_name;
247  mesh_name << "cond.mesh";
248  cond_name << "cond_j.gf";
249 
250  ofstream mesh_ofs(mesh_name.str().c_str());
251  mesh_ofs.precision(8);
252  mesh_cond.Print(mesh_ofs);
253 
254  ofstream cond_ofs(cond_name.str().c_str());
255  cond_ofs.precision(8);
256  j_cond.Save(cond_ofs);
257  }
258  // 6b. Send the current density, computed on the SubMesh, to a GLVis server.
259  if (visualization)
260  {
261  char vishost[] = "localhost";
262  int visport = 19916;
263  socketstream port_sock(vishost, visport);
264  port_sock.precision(8);
265  port_sock << "solution\n" << mesh_cond << j_cond
266  << "window_title 'Conductor J'"
267  << "window_geometry 400 0 400 350" << flush;
268  }
269 
270  // 7. Define a parallel finite element space on the full mesh. Here we use
271  // the H(curl) finite elements for the vector potential and H(div) for the
272  // current density.
273  ND_FECollection fec_nd(order, dim);
274  RT_FECollection fec_rt(order - 1, dim);
275  FiniteElementSpace fespace_nd(&mesh, &fec_nd);
276  FiniteElementSpace fespace_rt(&mesh, &fec_rt);
277 
278  GridFunction j_full(&fespace_rt);
279  j_full = 0.0;
280  mesh_cond.Transfer(j_cond, j_full);
281 
282  // 7a. Send the transferred current density to a GLVis server.
283  if (visualization)
284  {
285  char vishost[] = "localhost";
286  int visport = 19916;
287  socketstream sol_sock(vishost, visport);
288  sol_sock.precision(8);
289  sol_sock << "solution\n" << mesh << j_full
290  << "window_title 'J Full'"
291  << "window_geometry 400 430 400 350" << flush;
292  }
293 
294  // 8. Determine the list of true (i.e. parallel conforming) essential
295  // boundary dofs. In this example, the boundary conditions are defined by
296  // marking all the boundary attributes except for those on a symmetry
297  // plane as essential (Dirichlet) and converting them to a list of true
298  // dofs.
299  Array<int> ess_tdof_list;
300  Array<int> ess_bdr;
301  if (mesh.bdr_attributes.Size())
302  {
303  ess_bdr.SetSize(mesh.bdr_attributes.Max());
304  ess_bdr = 1;
305  for (int i=0; i<sym_plane_attr.Size(); i++)
306  {
307  ess_bdr[sym_plane_attr[i]-1] = 0;
308  }
309  fespace_nd.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
310  }
311 
312  // 9. Set up the parallel linear form b(.) which corresponds to the
313  // right-hand side of the FEM linear system, which in this case is
314  // (J,W_i) where J is given by the function H(div) field transferred
315  // from the SubMesh and W_i are the basis functions in the finite
316  // element fespace.
317  VectorGridFunctionCoefficient jCoef(&j_full);
318  LinearForm b(&fespace_nd);
319  b.AddDomainIntegrator(new VectorFEDomainLFIntegrator(jCoef));
320  b.Assemble();
321 
322  // 10. Define the solution vector x as a parallel finite element grid
323  // function corresponding to fespace. Initialize x to zero.
324  GridFunction x(&fespace_nd);
325  x = 0.0;
326 
327  // 11. Set up the parallel bilinear form corresponding to the EM diffusion
328  // operator curl muinv curl + delta I, by adding the curl-curl and the
329  // mass domain integrators. For standard magnetostatics equations choose
330  // delta << 1. Larger values of delta should make the linear system
331  // easier to solve at the expense of resembling a diffusive quasistatic
332  // magnetic field. A reasonable balance must be found whenever the mesh
333  // or problem setup is altered.
334  ConstantCoefficient muinv(1.0);
335  ConstantCoefficient delta(delta_const);
336  BilinearForm a(&fespace_nd);
337  if (pa_) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
338  a.AddDomainIntegrator(new CurlCurlIntegrator(muinv));
339  a.AddDomainIntegrator(new VectorFEMassIntegrator(delta));
340 
341  // 12. Assemble the parallel bilinear form and the corresponding linear
342  // system, applying any necessary transformations such as: parallel
343  // assembly, eliminating boundary conditions, applying conforming
344  // constraints for non-conforming AMR, static condensation, etc.
345  if (static_cond) { a.EnableStaticCondensation(); }
346  a.Assemble();
347 
348  OperatorPtr A;
349  Vector B, X;
350  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
351 
352  // 13. Solve the system AX=B
353  if (pa_) // Jacobi preconditioning in partial assembly mode
354  {
355  cout << "\nSolving for magnetic vector potential "
356  << "using CG with a Jacobi preconditioner" << endl;
357 
358  OperatorJacobiSmoother M(a, ess_tdof_list);
359  PCG(*A, M, B, X, 1, 1000, 1e-12, 0.0);
360  }
361  else
362  {
363 #ifndef MFEM_USE_SUITESPARSE
364  cout << "\nSolving for magnetic vector potential "
365  << "using CG with a Gauss-Seidel preconditioner" << endl;
366 
367  // 13a. Define a simple symmetric Gauss-Seidel preconditioner and use
368  // it to solve the system Ax=b with PCG.
369  GSSmoother M((SparseMatrix&)(*A));
370  PCG(*A, M, B, X, 1, 500, 1e-12, 0.0);
371 #else
372  cout << "\nSolving for magnetic vector potential "
373  << "using UMFPack" << endl;
374 
375  // 13a. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the
376  // system.
377  UMFPackSolver umf_solver;
378  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
379  umf_solver.SetOperator(*A);
380  umf_solver.Mult(B, X);
381 #endif
382  }
383 
384  // 14. Recover the parallel grid function corresponding to X. This is the
385  // local finite element solution on each processor.
386  a.RecoverFEMSolution(X, b, x);
387 
388  // 15. Save the refined mesh and the solution in parallel. This output can
389  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
390  {
391  ostringstream mesh_name, sol_name;
392  mesh_name << "refined.mesh";
393  sol_name << "sol.gf";
394 
395  ofstream mesh_ofs(mesh_name.str().c_str());
396  mesh_ofs.precision(8);
397  mesh.Print(mesh_ofs);
398 
399  ofstream sol_ofs(sol_name.str().c_str());
400  sol_ofs.precision(8);
401  x.Save(sol_ofs);
402  }
403 
404  // 16. Send the solution by socket to a GLVis server.
405  if (visualization)
406  {
407  char vishost[] = "localhost";
408  int visport = 19916;
409  socketstream sol_sock(vishost, visport);
410  sol_sock.precision(8);
411  sol_sock << "solution\n" << mesh << x
412  << "window_title 'Vector Potential'"
413  << "window_geometry 800 0 400 350" << flush;
414  }
415 
416  // 17. Compute the magnetic flux as the curl of the solution
417  DiscreteLinearOperator curl(&fespace_nd, &fespace_rt);
419  curl.Assemble();
420  curl.Finalize();
421 
422  GridFunction dx(&fespace_rt);
423  curl.Mult(x, dx);
424 
425  // 18. Save the curl of the solution in parallel. This output can be viewed
426  // later using GLVis: "glvis -np <np> -m mesh -g dsol".
427  {
428  ostringstream dsol_name;
429  dsol_name << "dsol.gf";
430 
431  ofstream dsol_ofs(dsol_name.str().c_str());
432  dsol_ofs.precision(8);
433  dx.Save(dsol_ofs);
434  }
435 
436  // 19. Send the curl of the solution by socket to a GLVis server.
437  if (visualization)
438  {
439  char vishost[] = "localhost";
440  int visport = 19916;
441  socketstream sol_sock(vishost, visport);
442  sol_sock.precision(8);
443  sol_sock << "solution\n" << mesh << dx
444  << "window_title 'Magnetic Flux'"
445  << "window_geometry 1200 0 400 350" << flush;
446  }
447 
448  // 20. Clean exit
449  return 0;
450 }
451 
453  const Array<int> &phi0_attr,
454  const Array<int> &phi1_attr,
455  const Array<int> &jn_zero_attr,
456  GridFunction &j_cond)
457 {
458  // Extract the finite element space and mesh on which j_cond is defined
459  FiniteElementSpace &fes_cond_rt = *j_cond.FESpace();
460  Mesh &mesh_cond = *fes_cond_rt.GetMesh();
461  int dim = mesh_cond.Dimension();
462 
463  // Define a parallel finite element space on the SubMesh. Here we use the H1
464  // finite elements for the electrostatic potential.
465  H1_FECollection fec_h1(order, dim);
466  FiniteElementSpace fes_cond_h1(&mesh_cond, &fec_h1);
467 
468  // Define the conductivity coefficient and the boundaries associated with the
469  // fixed potentials phi0 and phi1 which will drive the current.
470  ConstantCoefficient sigmaCoef(1.0);
471  Array<int> ess_bdr_phi(mesh_cond.bdr_attributes.Max());
472  Array<int> ess_bdr_j(mesh_cond.bdr_attributes.Max());
473  Array<int> ess_bdr_tdof_phi;
474  ess_bdr_phi = 0;
475  ess_bdr_j = 0;
476  for (int i=0; i<phi0_attr.Size(); i++)
477  {
478  ess_bdr_phi[phi0_attr[i]-1] = 1;
479  }
480  for (int i=0; i<phi1_attr.Size(); i++)
481  {
482  ess_bdr_phi[phi1_attr[i]-1] = 1;
483  }
484  for (int i=0; i<jn_zero_attr.Size(); i++)
485  {
486  ess_bdr_j[jn_zero_attr[i]-1] = 1;
487  }
488  fes_cond_h1.GetEssentialTrueDofs(ess_bdr_phi, ess_bdr_tdof_phi);
489 
490  // Setup the bilinear form corresponding to -Div(sigma Grad phi)
491  BilinearForm a_h1(&fes_cond_h1);
492  a_h1.AddDomainIntegrator(new DiffusionIntegrator(sigmaCoef));
493  a_h1.Assemble();
494 
495  // Set the r.h.s. to zero
496  LinearForm b_h1(&fes_cond_h1);
497  b_h1 = 0.0;
498 
499  // Setup the boundary conditions on phi
500  ConstantCoefficient one(1.0);
501  ConstantCoefficient zero(0.0);
502  GridFunction phi_h1(&fes_cond_h1);
503  phi_h1 = 0.0;
504 
505  Array<int> bdr0(mesh_cond.bdr_attributes.Max()); bdr0 = 0;
506  for (int i=0; i<phi0_attr.Size(); i++)
507  {
508  bdr0[phi0_attr[i]-1] = 1;
509  }
510  phi_h1.ProjectBdrCoefficient(zero, bdr0);
511 
512  Array<int> bdr1(mesh_cond.bdr_attributes.Max()); bdr1 = 0;
513  for (int i=0; i<phi1_attr.Size(); i++)
514  {
515  bdr1[phi1_attr[i]-1] = 1;
516  }
517  phi_h1.ProjectBdrCoefficient(one, bdr1);
518 
519  {
520  OperatorPtr A;
521  Vector B, X;
522  a_h1.FormLinearSystem(ess_bdr_tdof_phi, phi_h1, b_h1, A, X, B);
523 
524  // Solve the linear system
525  if (!pa_)
526  {
527 #ifndef MFEM_USE_SUITESPARSE
528  cout << "\nSolving for electric potential using PCG "
529  << "with a Gauss-Seidel preconditioner" << endl;
530 
531  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
532  GSSmoother M((SparseMatrix&)(*A));
533  PCG(*A, M, B, X, 1, 200, 1e-12, 0.0);
534 #else
535  cout << "\nSolving for electric potential using UMFPack" << endl;
536 
537  // If MFEM was compiled with SuiteSparse,
538  // use UMFPACK to solve the system.
539  UMFPackSolver umf_solver;
540  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
541  umf_solver.SetOperator(*A);
542  umf_solver.Mult(B, X);
543 #endif
544  }
545  else
546  {
547  cout << "\nSolving for electric potential using CG" << endl;
548 
549  if (UsesTensorBasis(fes_cond_h1))
550  {
551  if (algebraic_ceed_)
552  {
553  ceed::AlgebraicSolver M(a_h1, ess_bdr_tdof_phi);
554  PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
555  }
556  else
557  {
558  OperatorJacobiSmoother M(a_h1, ess_bdr_tdof_phi);
559  PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
560  }
561  }
562  else
563  {
564  CG(*A, B, X, 1, 400, 1e-12, 0.0);
565  }
566  }
567  a_h1.RecoverFEMSolution(X, b_h1, phi_h1);
568  }
569 
570  {
571  char vishost[] = "localhost";
572  int visport = 19916;
573  socketstream port_sock(vishost, visport);
574  port_sock.precision(8);
575  port_sock << "solution\n" << mesh_cond << phi_h1
576  << "window_title 'Conductor Potential'"
577  << "window_geometry 0 0 400 350" << flush;
578  }
579 
580  // Solve for the current density J = -sigma Grad phi with boundary conditions
581  // J.n = 0 on the walls of the conductor but not on the ports where phi=0 and
582  // phi=1.
583 
584  // J will be computed in H(div) so we need an RT mass matrix
585  BilinearForm m_rt(&fes_cond_rt);
587  m_rt.Assemble();
588 
589  // Assemble the (sigma Grad phi) operator
590  MixedBilinearForm d_h1(&fes_cond_h1, &fes_cond_rt);
592  d_h1.Assemble();
593 
594  // Compute the r.h.s, b_rt = sigma E = -sigma Grad phi
595  LinearForm b_rt(&fes_cond_rt);
596  d_h1.Mult(phi_h1, b_rt);
597  b_rt *= -1.0;
598 
599  // Apply the necessary boundary conditions and solve for J in H(div)
600  cout << "\nSolving for current density in H(Div) "
601  << "using diagonally scaled CG" << endl;
602  cout << "Size of linear system: "
603  << fes_cond_rt.GetTrueVSize() << endl;
604 
605  Array<int> ess_bdr_tdof_rt;
606  OperatorPtr M;
607  Vector B, X;
608 
609  fes_cond_rt.GetEssentialTrueDofs(ess_bdr_j, ess_bdr_tdof_rt);
610 
611  j_cond = 0.0;
612  m_rt.FormLinearSystem(ess_bdr_tdof_rt, j_cond, b_rt, M, X, B);
613 
614  CGSolver cg;
615  cg.SetRelTol(1e-12);
616  cg.SetMaxIter(2000);
617  cg.SetPrintLevel(1);
618  cg.SetOperator(*M);
619  cg.Mult(B, X);
620  m_rt.RecoverFEMSolution(X, b_rt, j_cond);
621 }
int visport
Conjugate gradient method.
Definition: solvers.hpp:493
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
Integrator for (curl u, curl v) for Nedelec elements.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:718
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:587
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
STL namespace.
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization.
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition: fespace.hpp:1306
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
Data type for Gauss-Seidel smoother of sparse matrix.
int main(int argc, char *argv[])
Definition: ex34.cpp:60
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:1070
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
Data type sparse matrix.
Definition: sparsemat.hpp:50
char vishost[]
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
Definition: mesh.cpp:1572
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:302
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
double delta
Definition: lissajous.cpp:43
void Assemble(int skip_zeros=1)
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:898
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:913
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:691
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:380
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:712
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:555
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void SetRelTol(double rtol)
Definition: solvers.hpp:199
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:1081
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
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:687
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3194
Wrapper for AlgebraicMultigrid object.
Definition: algebraic.hpp:185
double a
Definition: lissajous.cpp:41
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
int dim
Definition: ex24.cpp:53
static void Transfer(const GridFunction &src, GridFunction &dst)
Transfer the dofs of a GridFunction.
Definition: submesh.cpp:198
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Subdomain representation of a topological parent in another Mesh.
Definition: submesh.hpp:42
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition: lininteg.hpp:346
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition: mesh.hpp:1193
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:454
Vector data type.
Definition: vector.hpp:58
void ComputeCurrentDensityOnSubMesh(int order, const Array< int > &phi0_attr, const Array< int > &phi1_attr, const Array< int > &jn_zero_attr, GridFunction &j_cond)
Definition: ex34.cpp:452
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:2011
Vector coefficient defined by a vector GridFunction.
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:24
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3696
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition: mesh.hpp:273
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:468
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:3099