MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
ex34p.cpp
Go to the documentation of this file.
1// MFEM Example 34 - Parallel Version
2//
3// Compile with: make ex34p
4//
5// Sample runs: mpirun -np 4 ex34p -o 2
6// mpirun -np 4 ex34p -o 2 -hex -pa
7//
8// Device sample runs:
9// mpirun -np 4 ex34p -o 2 -hex -pa -d cuda
10// mpirun -np 4 ex34p -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
48using namespace std;
49using namespace mfem;
50
52 bool visualization,
53 const Array<int> &phi0_attr,
54 const Array<int> &phi1_attr,
55 const Array<int> &jn_zero_attr,
56 ParGridFunction &j_cond);
57
58int main(int argc, char *argv[])
59{
60 // 1. Initialize MPI and HYPRE.
61 Mpi::Init(argc, argv);
62 int num_procs = Mpi::WorldSize();
63 int myid = Mpi::WorldRank();
65
66 // 2. Parse command-line options.
67 const char *mesh_file = "../data/fichera-mixed.mesh";
68 Array<int> cond_attr;
69 Array<int> submesh_elems;
70 Array<int> sym_plane_attr;
71 Array<int> phi0_attr;
72 Array<int> phi1_attr;
73 Array<int> jn_zero_attr;
74 int ser_ref_levels = 1;
75 int par_ref_levels = 1;
76 int order = 1;
77 real_t delta_const = 1e-6;
78 bool mixed = true;
79 bool static_cond = false;
80 bool pa = false;
81 const char *device_config = "cpu";
82 bool visualization = true;
83#ifdef MFEM_USE_AMGX
84 bool useAmgX = false;
85#endif
86
87 OptionsParser args(argc, argv);
88 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
89 "Number of times to refine the mesh uniformly in serial.");
90 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
91 "Number of times to refine the mesh uniformly in parallel.");
92 args.AddOption(&order, "-o", "--order",
93 "Finite element order (polynomial degree).");
94 args.AddOption(&delta_const, "-mc", "--magnetic-cond",
95 "Magnetic Conductivity");
96 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
97 "--no-static-condensation", "Enable static condensation.");
98 args.AddOption(&mixed, "-mixed", "--mixed-mesh", "-hex",
99 "--hex-mesh", "Mixed mesh of hexahedral mesh.");
100 args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
101 "--no-partial-assembly", "Enable Partial Assembly.");
102 args.AddOption(&device_config, "-d", "--device",
103 "Device configuration string, see Device::Configure().");
104 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
105 "--no-visualization",
106 "Enable or disable GLVis visualization.");
107#ifdef MFEM_USE_AMGX
108 args.AddOption(&useAmgX, "-amgx", "--useAmgX", "-no-amgx",
109 "--no-useAmgX",
110 "Enable or disable AmgX in MatrixFreeAMS.");
111#endif
112
113 args.Parse();
114 if (!args.Good())
115 {
116 if (myid == 0)
117 {
118 args.PrintUsage(cout);
119 }
120 return 1;
121 }
122 if (myid == 0)
123 {
124 args.PrintOptions(cout);
125 }
126
127 if (!mixed || pa)
128 {
129 mesh_file = "../data/fichera.mesh";
130 }
131
132 if (submesh_elems.Size() == 0)
133 {
134 if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0)
135 {
136 submesh_elems.SetSize(5);
137 submesh_elems[0] = 0;
138 submesh_elems[1] = 2;
139 submesh_elems[2] = 3;
140 submesh_elems[3] = 4;
141 submesh_elems[4] = 9;
142 }
143 else if (strcmp(mesh_file, "../data/fichera.mesh") == 0)
144 {
145 submesh_elems.SetSize(7);
146 submesh_elems[0] = 10;
147 submesh_elems[1] = 14;
148 submesh_elems[2] = 34;
149 submesh_elems[3] = 36;
150 submesh_elems[4] = 37;
151 submesh_elems[5] = 38;
152 submesh_elems[6] = 39;
153 }
154 }
155 if (sym_plane_attr.Size() == 0)
156 {
157 if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
158 strcmp(mesh_file, "../data/fichera.mesh") == 0)
159 {
160 sym_plane_attr.SetSize(8);
161 sym_plane_attr[0] = 9;
162 sym_plane_attr[1] = 10;
163 sym_plane_attr[2] = 11;
164 sym_plane_attr[3] = 12;
165 sym_plane_attr[4] = 13;
166 sym_plane_attr[5] = 14;
167 sym_plane_attr[6] = 15;
168 sym_plane_attr[7] = 16;
169 }
170 }
171 if (phi0_attr.Size() == 0)
172 {
173 if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
174 strcmp(mesh_file, "../data/fichera.mesh") == 0)
175 {
176 phi0_attr.Append(2);
177 }
178 }
179 if (phi1_attr.Size() == 0)
180 {
181 if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
182 strcmp(mesh_file, "../data/fichera.mesh") == 0)
183 {
184 phi1_attr.Append(23);
185 }
186 }
187 if (jn_zero_attr.Size() == 0)
188 {
189 if (strcmp(mesh_file, "../data/fichera-mixed.mesh") == 0 ||
190 strcmp(mesh_file, "../data/fichera.mesh") == 0)
191 {
192 jn_zero_attr.Append(25);
193 }
194 for (int i=0; i<sym_plane_attr.Size(); i++)
195 {
196 jn_zero_attr.Append(sym_plane_attr[i]);
197 }
198 }
199
200 // 3. Enable hardware devices such as GPUs, and programming models such as
201 // CUDA, OCCA, RAJA and OpenMP based on command line options.
202 Device device(device_config);
203 if (myid == 0) { device.Print(); }
204
205 // 4. Read the (serial) mesh from the given mesh file on all processors. We
206 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
207 // and volume meshes with the same code.
208 Mesh *mesh = new Mesh(mesh_file, 1, 1);
209 int dim = mesh->Dimension();
210
211 if (!mixed || pa)
212 {
213 mesh->UniformRefinement();
214
215 if (ser_ref_levels > 0)
216 {
217 ser_ref_levels--;
218 }
219 else
220 {
221 par_ref_levels--;
222 }
223 }
224
225 int submesh_attr = -1;
226 if (cond_attr.Size() == 0 && submesh_elems.Size() > 0)
227 {
228 int max_attr = mesh->attributes.Max();
229 submesh_attr = max_attr + 1;
230
231 for (int i=0; i<submesh_elems.Size(); i++)
232 {
233 mesh->SetAttribute(submesh_elems[i], submesh_attr);
234 }
235 mesh->SetAttributes();
236
237 if (cond_attr.Size() == 0)
238 {
239 cond_attr.Append(submesh_attr);
240 }
241 }
242
243 // 5. Refine the serial mesh on all processors to increase the resolution. In
244 // this example we do 'ref_levels' of uniform refinement.
245 {
246 int ref_levels = ser_ref_levels;
247 for (int l = 0; l < ref_levels; l++)
248 {
249 mesh->UniformRefinement();
250 }
251 }
252
253 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
254 // this mesh further in parallel to increase the resolution. Once the
255 // parallel mesh is defined, the serial mesh can be deleted.
256 ParMesh pmesh(MPI_COMM_WORLD, *mesh);
257 delete mesh;
258 {
259 for (int l = 0; l < par_ref_levels; l++)
260 {
261 pmesh.UniformRefinement();
262 }
263 }
264
265 // 6b. Extract a submesh covering a portion of the domain
266 ParSubMesh pmesh_cond(ParSubMesh::CreateFromDomain(pmesh, cond_attr));
267
268 // 7. Define a suitable finite element space on the SubMesh and compute
269 // the current density as an H(div) field.
270 RT_FECollection fec_cond_rt(order - 1, dim);
271 ParFiniteElementSpace fes_cond_rt(&pmesh_cond, &fec_cond_rt);
272 ParGridFunction j_cond(&fes_cond_rt);
273
274 ComputeCurrentDensityOnSubMesh(order, visualization,
275 phi0_attr, phi1_attr, jn_zero_attr, j_cond);
276
277 // 7a. Save the SubMesh and associated current density in parallel. This
278 // output can be viewed later using GLVis:
279 // "glvis -np <np> -m cond_mesh -g cond_j"
280 {
281 ostringstream mesh_name, cond_name;
282 mesh_name << "cond_mesh." << setfill('0') << setw(6) << myid;
283 cond_name << "cond_j." << setfill('0') << setw(6) << myid;
284
285 ofstream mesh_ofs(mesh_name.str().c_str());
286 mesh_ofs.precision(8);
287 pmesh_cond.Print(mesh_ofs);
288
289 ofstream cond_ofs(cond_name.str().c_str());
290 cond_ofs.precision(8);
291 j_cond.Save(cond_ofs);
292 }
293
294 // 7b. Send the current density, computed on the SubMesh, to a GLVis server.
295 if (visualization)
296 {
297 char vishost[] = "localhost";
298 int visport = 19916;
299 socketstream port_sock(vishost, visport);
300 port_sock << "parallel " << num_procs << " " << myid << "\n";
301 port_sock.precision(8);
302 port_sock << "solution\n" << pmesh_cond << j_cond
303 << "window_title 'Conductor J'"
304 << "window_geometry 400 0 400 350" << flush;
305 }
306
307 // 8. Define a parallel finite element space on the full mesh. Here we use
308 // the H(curl) finite elements for the vector potential and H(div) for the
309 // current density.
310 ND_FECollection fec_nd(order, dim);
311 RT_FECollection fec_rt(order - 1, dim);
312 ParFiniteElementSpace fespace_nd(&pmesh, &fec_nd);
313 ParFiniteElementSpace fespace_rt(&pmesh, &fec_rt);
314
315 ParGridFunction j_full(&fespace_rt);
316 j_full = 0.0;
317 pmesh_cond.Transfer(j_cond, j_full);
318
319 // 8a. Send the transferred current density to a GLVis server.
320 if (visualization)
321 {
322 char vishost[] = "localhost";
323 int visport = 19916;
324 socketstream sol_sock(vishost, visport);
325 sol_sock << "parallel " << num_procs << " " << myid << "\n";
326 sol_sock.precision(8);
327 sol_sock << "solution\n" << pmesh << j_full
328 << "window_title 'J Full'"
329 << "window_geometry 400 430 400 350" << flush;
330 }
331
332 // 9. Determine the list of true (i.e. parallel conforming) essential
333 // boundary dofs. In this example, the boundary conditions are defined
334 // by marking all the boundary attributes except for those on a symmetry
335 // plane as essential (Dirichlet) and converting them to a list of
336 // true dofs.
337 Array<int> ess_tdof_list;
338 Array<int> ess_bdr;
339 if (pmesh.bdr_attributes.Size())
340 {
341 ess_bdr.SetSize(pmesh.bdr_attributes.Max());
342 ess_bdr = 1;
343 for (int i=0; i<sym_plane_attr.Size(); i++)
344 {
345 ess_bdr[sym_plane_attr[i]-1] = 0;
346 }
347 fespace_nd.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
348 }
349
350 // 10. Set up the parallel linear form b(.) which corresponds to the
351 // right-hand side of the FEM linear system, which in this case is
352 // (J,W_i) where J is given by the function H(div) field transferred
353 // from the SubMesh and W_i are the basis functions in the finite
354 // element fespace.
355 VectorGridFunctionCoefficient jCoef(&j_full);
356 ParLinearForm b(&fespace_nd);
357 b.AddDomainIntegrator(new VectorFEDomainLFIntegrator(jCoef));
358 b.Assemble();
359
360 // 11. Define the solution vector x as a parallel finite element grid
361 // function corresponding to fespace. Initialize x to zero.
362 ParGridFunction x(&fespace_nd);
363 x = 0.0;
364
365 // 12. Set up the parallel bilinear form corresponding to the EM diffusion
366 // operator curl muinv curl + delta I, by adding the curl-curl and the
367 // mass domain integrators. For standard magnetostatics equations choose
368 // delta << 1. Larger values of delta should make the linear system
369 // easier to solve at the expense of resembling a diffusive quasistatic
370 // magnetic field. A reasonable balance must be found whenever the mesh
371 // or problem setup is altered.
372 ConstantCoefficient muinv(1.0);
373 ConstantCoefficient delta(delta_const);
374 ParBilinearForm a(&fespace_nd);
375 if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
376 a.AddDomainIntegrator(new CurlCurlIntegrator(muinv));
377 a.AddDomainIntegrator(new VectorFEMassIntegrator(delta));
378
379 // 13. Assemble the parallel bilinear form and the corresponding linear
380 // system, applying any necessary transformations such as: parallel
381 // assembly, eliminating boundary conditions, applying conforming
382 // constraints for non-conforming AMR, static condensation, etc.
383 if (static_cond) { a.EnableStaticCondensation(); }
384 a.Assemble();
385
386 OperatorPtr A;
387 Vector B, X;
388 a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
389
390 if (myid == 0)
391 {
392 cout << "\nSolving for magnetic vector potential "
393 << "using CG with AMS" << endl;
394 }
395
396 // 14. Solve the system AX=B using PCG with an AMS preconditioner.
397 if (pa)
398 {
399#ifdef MFEM_USE_AMGX
400 MatrixFreeAMS ams(a, *A, fespace_nd, &muinv, &delta, NULL, ess_bdr,
401 useAmgX);
402#else
403 MatrixFreeAMS ams(a, *A, fespace_nd, &muinv, &delta, NULL, ess_bdr);
404#endif
405 CGSolver cg(MPI_COMM_WORLD);
406 cg.SetRelTol(1e-12);
407 cg.SetMaxIter(1000);
408 cg.SetPrintLevel(1);
409 cg.SetOperator(*A);
410 cg.SetPreconditioner(ams);
411 cg.Mult(B, X);
412 }
413 else
414 {
415 if (myid == 0)
416 {
417 cout << "Size of linear system: "
418 << A.As<HypreParMatrix>()->GetGlobalNumRows() << endl;
419 }
420
421 ParFiniteElementSpace *prec_fespace =
422 (a.StaticCondensationIsEnabled() ? a.SCParFESpace() : &fespace_nd);
423 HypreAMS ams(*A.As<HypreParMatrix>(), prec_fespace);
424 HyprePCG pcg(*A.As<HypreParMatrix>());
425 pcg.SetTol(1e-12);
426 pcg.SetMaxIter(500);
427 pcg.SetPrintLevel(2);
428 pcg.SetPreconditioner(ams);
429 pcg.Mult(B, X);
430 }
431
432 // 15. Recover the parallel grid function corresponding to X. This is the
433 // local finite element solution on each processor.
434 a.RecoverFEMSolution(X, b, x);
435
436 // 16. Save the refined mesh and the solution in parallel. This output can
437 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
438 {
439 ostringstream mesh_name, sol_name;
440 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
441 sol_name << "sol." << setfill('0') << setw(6) << myid;
442
443 ofstream mesh_ofs(mesh_name.str().c_str());
444 mesh_ofs.precision(8);
445 pmesh.Print(mesh_ofs);
446
447 ofstream sol_ofs(sol_name.str().c_str());
448 sol_ofs.precision(8);
449 x.Save(sol_ofs);
450 }
451
452 // 17. Send the solution by socket to a GLVis server.
453 if (visualization)
454 {
455 char vishost[] = "localhost";
456 int visport = 19916;
457 socketstream sol_sock(vishost, visport);
458 sol_sock << "parallel " << num_procs << " " << myid << "\n";
459 sol_sock.precision(8);
460 sol_sock << "solution\n" << pmesh << x
461 << "window_title 'Vector Potential'"
462 << "window_geometry 800 0 400 350" << flush;
463 }
464
465 // 18. Compute the magnetic flux as the curl of the solution
466 ParDiscreteLinearOperator curl(&fespace_nd, &fespace_rt);
468 curl.Assemble();
469 curl.Finalize();
470
471 ParGridFunction dx(&fespace_rt);
472 curl.Mult(x, dx);
473
474 // 19. Save the curl of the solution in parallel. This output can
475 // be viewed later using GLVis: "glvis -np <np> -m mesh -g dsol".
476 {
477 ostringstream dsol_name;
478 dsol_name << "dsol." << setfill('0') << setw(6) << myid;
479
480 ofstream dsol_ofs(dsol_name.str().c_str());
481 dsol_ofs.precision(8);
482 dx.Save(dsol_ofs);
483 }
484
485 // 20. Send the curl of the solution by socket to a GLVis server.
486 if (visualization)
487 {
488 char vishost[] = "localhost";
489 int visport = 19916;
490 socketstream sol_sock(vishost, visport);
491 sol_sock << "parallel " << num_procs << " " << myid << "\n";
492 sol_sock.precision(8);
493 sol_sock << "solution\n" << pmesh << dx
494 << "window_title 'Magnetic Flux'"
495 << "window_geometry 1200 0 400 350" << flush;
496 }
497
498 // 21. Clean exit
499 return 0;
500}
501
503 bool visualization,
504 const Array<int> &phi0_attr,
505 const Array<int> &phi1_attr,
506 const Array<int> &jn_zero_attr,
507 ParGridFunction &j_cond)
508{
509 // Extract the finite element space and mesh on which j_cond is defined
510 ParFiniteElementSpace &fes_cond_rt = *j_cond.ParFESpace();
511 ParMesh &pmesh_cond = *fes_cond_rt.GetParMesh();
512 int myid = fes_cond_rt.GetMyRank();
513 int dim = pmesh_cond.Dimension();
514
515 // Define a parallel finite element space on the SubMesh. Here we use the
516 // H1 finite elements for the electrostatic potential.
517 H1_FECollection fec_h1(order, dim);
518 ParFiniteElementSpace fes_cond_h1(&pmesh_cond, &fec_h1);
519
520 // Define the conductivity coefficient and the boundaries associated with the
521 // fixed potentials phi0 and phi1 which will drive the current.
522 ConstantCoefficient sigmaCoef(1.0);
523 Array<int> ess_bdr_phi(pmesh_cond.bdr_attributes.Max());
524 Array<int> ess_bdr_j(pmesh_cond.bdr_attributes.Max());
525 Array<int> ess_bdr_tdof_phi;
526 ess_bdr_phi = 0;
527 ess_bdr_j = 0;
528 for (int i=0; i<phi0_attr.Size(); i++)
529 {
530 ess_bdr_phi[phi0_attr[i]-1] = 1;
531 }
532 for (int i=0; i<phi1_attr.Size(); i++)
533 {
534 ess_bdr_phi[phi1_attr[i]-1] = 1;
535 }
536 for (int i=0; i<jn_zero_attr.Size(); i++)
537 {
538 ess_bdr_j[jn_zero_attr[i]-1] = 1;
539 }
540 fes_cond_h1.GetEssentialTrueDofs(ess_bdr_phi, ess_bdr_tdof_phi);
541
542 // Setup the bilinear form corresponding to -Div(sigma Grad phi)
543 ParBilinearForm a_h1(&fes_cond_h1);
544 a_h1.AddDomainIntegrator(new DiffusionIntegrator(sigmaCoef));
545 a_h1.Assemble();
546
547 // Set the r.h.s. to zero
548 ParLinearForm b_h1(&fes_cond_h1);
549 b_h1 = 0.0;
550
551 // Setup the boundary conditions on phi
552 ConstantCoefficient one(1.0);
553 ConstantCoefficient zero(0.0);
554 ParGridFunction phi_h1(&fes_cond_h1);
555 phi_h1 = 0.0;
556
557 Array<int> bdr0(pmesh_cond.bdr_attributes.Max()); bdr0 = 0;
558 for (int i=0; i<phi0_attr.Size(); i++)
559 {
560 bdr0[phi0_attr[i]-1] = 1;
561 }
562 phi_h1.ProjectBdrCoefficient(zero, bdr0);
563
564 Array<int> bdr1(pmesh_cond.bdr_attributes.Max()); bdr1 = 0;
565 for (int i=0; i<phi1_attr.Size(); i++)
566 {
567 bdr1[phi1_attr[i]-1] = 1;
568 }
569 phi_h1.ProjectBdrCoefficient(one, bdr1);
570
571 // Solve the linear system using algebraic multigrid
572 {
573 if (myid == 0)
574 {
575 cout << "\nSolving for electric potential "
576 << "using CG with AMG" << endl;
577 }
578 OperatorPtr A;
579 Vector B, X;
580 a_h1.FormLinearSystem(ess_bdr_tdof_phi, phi_h1, b_h1, A, X, B);
581
582 HypreBoomerAMG prec;
583 CGSolver cg(MPI_COMM_WORLD);
584 cg.SetRelTol(1e-12);
585 cg.SetMaxIter(2000);
586 cg.SetPrintLevel(1);
587 cg.SetPreconditioner(prec);
588 cg.SetOperator(*A);
589 cg.Mult(B, X);
590 a_h1.RecoverFEMSolution(X, b_h1, phi_h1);
591 }
592
593 if (visualization)
594 {
595 int num_procs = fes_cond_h1.GetNRanks();
596 char vishost[] = "localhost";
597 int visport = 19916;
598 socketstream port_sock(vishost, visport);
599 port_sock << "parallel " << num_procs << " " << myid << "\n";
600 port_sock.precision(8);
601 port_sock << "solution\n" << pmesh_cond << phi_h1
602 << "window_title 'Conductor Potential'"
603 << "window_geometry 0 0 400 350" << flush;
604 }
605
606 // Solve for the current density J = -sigma Grad phi with boundary conditions
607 // J.n = 0 on the walls of the conductor but not on the ports where phi=0 and
608 // phi=1.
609
610 // J will be computed in H(div) so we need an RT mass matrix
611 ParBilinearForm m_rt(&fes_cond_rt);
613 m_rt.Assemble();
614
615 // Assemble the (sigma Grad phi) operator
616 ParMixedBilinearForm d_h1(&fes_cond_h1, &fes_cond_rt);
618 d_h1.Assemble();
619
620 // Compute the r.h.s, b_rt = sigma E = -sigma Grad phi
621 ParLinearForm b_rt(&fes_cond_rt);
622 d_h1.Mult(phi_h1, b_rt);
623 b_rt *= -1.0;
624
625 // Apply the necessary boundary conditions and solve for J in H(div)
626 HYPRE_BigInt glb_size_rt = fes_cond_rt.GlobalTrueVSize();
627 if (myid == 0)
628 {
629 cout << "\nSolving for current density in H(Div) "
630 << "using diagonally scaled CG" << endl;
631 cout << "Size of linear system: "
632 << glb_size_rt << endl;
633 }
634 Array<int> ess_bdr_tdof_rt;
635 OperatorPtr M;
636 Vector B, X;
637
638 fes_cond_rt.GetEssentialTrueDofs(ess_bdr_j, ess_bdr_tdof_rt);
639
640 j_cond = 0.0;
641 m_rt.FormLinearSystem(ess_bdr_tdof_rt, j_cond, b_rt, M, X, B);
642
643 HypreDiagScale prec;
644
645 CGSolver cg(MPI_COMM_WORLD);
646 cg.SetRelTol(1e-12);
647 cg.SetMaxIter(2000);
648 cg.SetPrintLevel(1);
649 cg.SetPreconditioner(prec);
650 cg.SetOperator(*M);
651 cg.Mult(B, X);
652 m_rt.RecoverFEMSolution(X, b_rt, j_cond);
653}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
A coefficient that is constant across space and time.
Integrator for for Nedelec elements.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
virtual void Assemble(int skip_zeros=1)
Construct the internal matrix representation of the discrete linear operator.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The Auxiliary-space Maxwell Solver in hypre.
Definition hypre.hpp:1845
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
Jacobi preconditioner in hypre.
Definition hypre.hpp:1481
PCG solver in hypre.
Definition hypre.hpp:1275
void SetPrintLevel(int print_lvl)
Definition hypre.cpp:4156
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition hypre.cpp:4161
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre's PCG.
Definition hypre.cpp:4184
void SetMaxIter(int max_iter)
Definition hypre.cpp:4146
void SetTol(real_t tol)
Definition hypre.cpp:4136
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
An auxiliary Maxwell solver for a high-order curl-curl system without high-order assembly.
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition mesh.hpp:1336
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition mesh.hpp:280
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
Definition mesh.cpp:1604
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY.
void Assemble(int skip_zeros=1)
virtual void Mult(const Vector &x, Vector &y) const
Matrix multiplication: .
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
static int WorldRank()
Return the MPI rank in MPI_COMM_WORLD.
static int WorldSize()
Return the size of MPI_COMM_WORLD.
static void Init(int &argc, char **&argv, int required=default_thread_required, int *provided=nullptr)
Singleton creation with Mpi::Init(argc, argv).
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition fe_coll.hpp:465
Pointer to an Operator of a specified type.
Definition handle.hpp:34
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Class for parallel bilinear form.
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(....
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalTrueVSize() const
Definition pfespace.hpp:285
ParMesh * GetParMesh() const
Definition pfespace.hpp:277
Class for parallel grid function.
Definition pgridfunc.hpp:33
void Save(std::ostream &out) const override
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
ParFiniteElementSpace * ParFESpace() const
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
Class for parallel bilinear form using different test and trial FE spaces.
Subdomain representation of a topological parent in another ParMesh.
Definition psubmesh.hpp:52
static ParSubMesh CreateFromDomain(const ParMesh &parent, Array< int > &domain_attributes)
Create a domain ParSubMesh from it's parent.
Definition psubmesh.cpp:26
static void Transfer(const ParGridFunction &src, ParGridFunction &dst)
Transfer the dofs of a ParGridFunction.
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
for VectorFiniteElements (Nedelec, Raviart-Thomas)
Definition lininteg.hpp:347
Vector coefficient defined by a vector GridFunction.
Vector data type.
Definition vector.hpp:80
int dim
Definition ex24.cpp:53
void ComputeCurrentDensityOnSubMesh(int order, bool visualization, const Array< int > &phi0_attr, const Array< int > &phi1_attr, const Array< int > &jn_zero_attr, ParGridFunction &j_cond)
Definition ex34p.cpp:502
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t delta
Definition lissajous.cpp:43
real_t a
Definition lissajous.cpp:41
const int visport
float real_t
Definition config.hpp:43
const char vishost[]