MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
ex1p.cpp
Go to the documentation of this file.
1// MFEM Example 1 - Parallel High-Performance Version
2//
3// Compile with: make ex1p
4//
5// Sample runs: mpirun -np 4 ex1p -m ../../data/fichera.mesh -perf -mf -pc lor
6// mpirun -np 4 ex1p -m ../../data/fichera.mesh -perf -asm -pc ho
7// mpirun -np 4 ex1p -m ../../data/fichera.mesh -perf -asm -pc ho -sc
8// mpirun -np 4 ex1p -m ../../data/fichera.mesh -std -asm -pc ho
9// mpirun -np 4 ex1p -m ../../data/fichera.mesh -std -asm -pc ho -sc
10// mpirun -np 4 ex1p -m ../../data/amr-hex.mesh -perf -asm -pc ho -sc
11// mpirun -np 4 ex1p -m ../../data/amr-hex.mesh -std -asm -pc ho -sc
12// mpirun -np 4 ex1p -m ../../data/ball-nurbs.mesh -perf -asm -pc ho -sc
13// mpirun -np 4 ex1p -m ../../data/ball-nurbs.mesh -std -asm -pc ho -sc
14// mpirun -np 4 ex1p -m ../../data/pipe-nurbs.mesh -perf -mf -pc lor
15// mpirun -np 4 ex1p -m ../../data/pipe-nurbs.mesh -std -asm -pc ho -sc
16// mpirun -np 4 ex1p -m ../../data/star.mesh -perf -mf -pc lor
17// mpirun -np 4 ex1p -m ../../data/star.mesh -perf -asm -pc ho
18// mpirun -np 4 ex1p -m ../../data/star.mesh -perf -asm -pc ho -sc
19// mpirun -np 4 ex1p -m ../../data/star.mesh -std -asm -pc ho
20// mpirun -np 4 ex1p -m ../../data/star.mesh -std -asm -pc ho -sc
21// mpirun -np 4 ex1p -m ../../data/amr-quad.mesh -perf -asm -pc ho -sc
22// mpirun -np 4 ex1p -m ../../data/amr-quad.mesh -std -asm -pc ho -sc
23// mpirun -np 4 ex1p -m ../../data/disc-nurbs.mesh -perf -asm -pc ho -sc
24// mpirun -np 4 ex1p -m ../../data/disc-nurbs.mesh -std -asm -pc ho -sc
25//
26// Description: This example code demonstrates the use of MFEM to define a
27// simple finite element discretization of the Poisson problem
28// -Delta u = 1 with homogeneous Dirichlet boundary conditions.
29// Specifically, we discretize using a FE space of the specified
30// order, or if order < 1 using an isoparametric/isogeometric
31// space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
32// NURBS mesh, etc.)
33//
34// The example highlights the use of mesh refinement, finite
35// element grid functions, as well as linear and bilinear forms
36// corresponding to the left-hand side and right-hand side of the
37// discrete linear system. We also cover the explicit elimination
38// of essential boundary conditions, static condensation, and the
39// optional connection to the GLVis tool for visualization.
40
41#include "mfem-performance.hpp"
42#include <fstream>
43#include <iostream>
44
45using namespace std;
46using namespace mfem;
47
48enum class PCType { NONE, LOR, HO };
49
50// Define template parameters for optimized build.
51template <int dim> struct geom_t { };
52template <>
53struct geom_t<2> { static const Geometry::Type value = Geometry::SQUARE; };
54template <>
55struct geom_t<3> { static const Geometry::Type value = Geometry::CUBE; };
56
57const int mesh_p = 3; // mesh curvature (default: 3)
58const int sol_p = 3; // solution order (default: 3)
59
60template <int dim>
61struct ex1_t
62{
63 static const Geometry::Type geom = geom_t<dim>::value;
64 static const int rdim = Geometry::Constants<geom>::Dimension;
65 static const int ir_order = 2*sol_p+rdim-1;
66
67 // Static mesh type
68 using mesh_fe_t = H1_FiniteElement<geom,mesh_p>;
69 using mesh_fes_t = H1_FiniteElementSpace<mesh_fe_t>;
70 using mesh_t = TMesh<mesh_fes_t>;
71
72 // Static solution finite element space type
73 using sol_fe_t = H1_FiniteElement<geom,sol_p>;
74 using sol_fes_t = H1_FiniteElementSpace<sol_fe_t>;
75
76 // Static quadrature, coefficient and integrator types
77 using int_rule_t = TIntegrationRule<geom,ir_order>;
78 using coeff_t = TConstantCoefficient<>;
80
82
83 static int run(Mesh *mesh, int ser_ref_levels, int par_ref_levels, int order,
84 int basis, bool static_cond, PCType pc_choice, bool perf,
85 bool matrix_free, bool visualization, int visport = 19916);
86};
87
88int main(int argc, char *argv[])
89{
90 // 1. Initialize MPI and HYPRE.
91 Mpi::Init(argc, argv);
92 int myid = Mpi::WorldRank();
94
95 // 2. Parse command-line options.
96#ifdef MFEM_HPC_EX1_2D
97 const char *mesh_file = "../../data/star.mesh";
98#else
99 const char *mesh_file = "../../data/fichera.mesh";
100#endif
101 int ser_ref_levels = -1;
102 int par_ref_levels = 1;
103 int order = sol_p;
104 const char *basis_type = "G"; // Gauss-Lobatto
105 bool static_cond = false;
106 const char *pc = "lor";
107 bool perf = true;
108 bool matrix_free = true;
109 int visport = 19916;
110 bool visualization = 1;
111
112 OptionsParser args(argc, argv);
113 args.AddOption(&mesh_file, "-m", "--mesh",
114 "Mesh file to use.");
115 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
116 "Number of times to refine the mesh uniformly in serial;"
117 " -1 = auto: <= 10,000 elements.");
118 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
119 "Number of times to refine the mesh uniformly in parallel.");
120 args.AddOption(&order, "-o", "--order",
121 "Finite element order (polynomial degree) or -1 for"
122 " isoparametric space.");
123 args.AddOption(&basis_type, "-b", "--basis-type",
124 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
125 args.AddOption(&perf, "-perf", "--hpc-version", "-std", "--standard-version",
126 "Enable high-performance, tensor-based, assembly/evaluation.");
127 args.AddOption(&matrix_free, "-mf", "--matrix-free", "-asm", "--assembly",
128 "Use matrix-free evaluation or efficient matrix assembly in "
129 "the high-performance version.");
130 args.AddOption(&pc, "-pc", "--preconditioner",
131 "Preconditioner: lor - low-order-refined (matrix-free) AMG, "
132 "ho - high-order (assembled) AMG, none.");
133 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
134 "--no-static-condensation", "Enable static condensation.");
135 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
136 "--no-visualization",
137 "Enable or disable GLVis visualization.");
138 args.AddOption(&visport, "-p", "--send-port", "Socket for GLVis.");
139 args.Parse();
140 if (!args.Good())
141 {
142 if (myid == 0)
143 {
144 args.PrintUsage(cout);
145 }
146 return 1;
147 }
148 if (static_cond && perf && matrix_free)
149 {
150 if (myid == 0)
151 {
152 cout << "\nStatic condensation can not be used with matrix-free"
153 " evaluation!\n" << endl;
154 }
155 return 2;
156 }
157 MFEM_VERIFY(perf || !matrix_free,
158 "--standard-version is not compatible with --matrix-free");
159 if (myid == 0)
160 {
161 args.PrintOptions(cout);
162 }
163
164 PCType pc_choice;
165 if (!strcmp(pc, "ho")) { pc_choice = PCType::HO; }
166 else if (!strcmp(pc, "lor")) { pc_choice = PCType::LOR; }
167 else if (!strcmp(pc, "none")) { pc_choice = PCType::NONE; }
168 else
169 {
170 mfem_error("Invalid Preconditioner specified");
171 return 3;
172 }
173
174 if (myid == 0)
175 {
176 cout << "\nMFEM SIMD width: " << MFEM_SIMD_BYTES/sizeof(double)
177 << " doubles\n" << endl;
178 }
179
180 // See class BasisType in fem/fe_coll.hpp for available basis types
181 int basis = BasisType::GetType(basis_type[0]);
182 if (myid == 0)
183 {
184 cout << "Using " << BasisType::Name(basis) << " basis ..." << endl;
185 }
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 = new Mesh(mesh_file, 1, 1);
191 int dim = mesh->Dimension();
192
193 if (dim == 2)
194 {
195 return ex1_t<2>::run(mesh, ser_ref_levels, par_ref_levels, order, basis,
196 static_cond, pc_choice, perf, matrix_free,
197 visualization, visport);
198 }
199 else if (dim == 3)
200 {
201 return ex1_t<3>::run(mesh, ser_ref_levels, par_ref_levels, order,
202 basis, static_cond, pc_choice, perf, matrix_free,
203 visualization, visport);
204 }
205 else
206 {
207 MFEM_ABORT("Dimension must be 2 or 3.")
208 }
209
210 return 0;
211}
212
213template <int dim>
214int ex1_t<dim>::run(Mesh *mesh, int ser_ref_levels, int par_ref_levels,
215 int order, int basis, bool static_cond, PCType pc_choice,
216 bool perf, bool matrix_free, bool visualization, int visport)
217{
218 int num_procs, myid;
219 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
220 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
221
222 // 4. Check if the optimized version matches the given mesh
223 if (perf)
224 {
225 if (myid == 0)
226 {
227 cout << "High-performance version using integration rule with "
228 << int_rule_t::qpts << " points ..." << endl;
229 }
230 if (!mesh_t::MatchesGeometry(*mesh))
231 {
232 if (myid == 0)
233 {
234 cout << "The given mesh does not match the optimized 'geom' parameter.\n"
235 << "Recompile with suitable 'geom' value." << endl;
236 }
237 delete mesh;
238 return 4;
239 }
240 else if (!mesh_t::MatchesNodes(*mesh))
241 {
242 if (myid == 0)
243 {
244 cout << "Switching the mesh curvature to match the "
245 << "optimized value (order " << mesh_p << ") ..." << endl;
246 }
247 mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
248 }
249 }
250
251 // 5. Refine the serial mesh on all processors to increase the resolution. In
252 // this example we do 'ref_levels' of uniform refinement. We choose
253 // 'ref_levels' to be the largest number that gives a final mesh with no
254 // more than 10,000 elements, or as specified on the command line with the
255 // option '--refine-serial'.
256 {
257 int ref_levels = (ser_ref_levels != -1) ? ser_ref_levels :
258 (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
259 for (int l = 0; l < ref_levels; l++)
260 {
261 mesh->UniformRefinement();
262 }
263 }
264 if (!perf && mesh->NURBSext)
265 {
266 const int new_mesh_p = std::min(sol_p, mesh_p);
267 if (myid == 0)
268 {
269 cout << "NURBS mesh: switching the mesh curvature to be "
270 << "min(sol_p, mesh_p) = " << new_mesh_p << " ..." << endl;
271 }
272 mesh->SetCurvature(new_mesh_p, false, -1, Ordering::byNODES);
273 }
274
275 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
276 // this mesh further in parallel to increase the resolution. Once the
277 // parallel mesh is defined, the serial mesh can be deleted.
278 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
279 delete mesh;
280 {
281 for (int l = 0; l < par_ref_levels; l++)
282 {
283 pmesh->UniformRefinement();
284 }
285 }
286 if (pmesh->MeshGenerator() & 1) // simplex mesh
287 {
288 MFEM_VERIFY(pc_choice != PCType::LOR, "triangle and tet meshes do not "
289 "support the LOR preconditioner yet");
290 }
291
292 // 7. Define a parallel finite element space on the parallel mesh. Here we
293 // use continuous Lagrange finite elements of the specified order. If
294 // order < 1, we instead use an isoparametric/isogeometric space.
296 if (order > 0)
297 {
298 fec = new H1_FECollection(order, dim, basis);
299 }
300 else if (pmesh->GetNodes())
301 {
302 fec = pmesh->GetNodes()->OwnFEC();
303 if (myid == 0)
304 {
305 cout << "Using isoparametric FEs: " << fec->Name() << endl;
306 }
307 }
308 else
309 {
310 fec = new H1_FECollection(order = 1, dim, basis);
311 }
312 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
313 HYPRE_BigInt size = fespace->GlobalTrueVSize();
314 if (myid == 0)
315 {
316 cout << "Number of finite element unknowns: " << size << endl;
317 }
318
319 ParMesh pmesh_lor;
320 FiniteElementCollection *fec_lor = NULL;
321 ParFiniteElementSpace *fespace_lor = NULL;
322 if (pc_choice == PCType::LOR)
323 {
324 int basis_lor = basis;
325 if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
326 pmesh_lor = ParMesh::MakeRefined(*pmesh, order, basis_lor);
327 fec_lor = new H1_FECollection(1, dim);
328 fespace_lor = new ParFiniteElementSpace(&pmesh_lor, fec_lor);
329 }
330
331 // 8. Check if the optimized version matches the given space
332 if (perf && !sol_fes_t::Matches(*fespace))
333 {
334 if (myid == 0)
335 {
336 cout << "The given order does not match the optimized parameter.\n"
337 << "Recompile with suitable 'sol_p' value." << endl;
338 }
339 delete fespace;
340 delete fec;
341 delete mesh;
342 return 5;
343 }
344
345 // 9. Determine the list of true (i.e. parallel conforming) essential
346 // boundary dofs. In this example, the boundary conditions are defined
347 // by marking all the boundary attributes from the mesh as essential
348 // (Dirichlet) and converting them to a list of true dofs.
349 Array<int> ess_tdof_list;
350 if (pmesh->bdr_attributes.Size())
351 {
352 Array<int> ess_bdr(pmesh->bdr_attributes.Max());
353 ess_bdr = 1;
354 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
355 }
356
357 // 10. Set up the parallel linear form b(.) which corresponds to the
358 // right-hand side of the FEM linear system, which in this case is
359 // (1,phi_i) where phi_i are the basis functions in fespace.
360 ParLinearForm *b = new ParLinearForm(fespace);
361 ConstantCoefficient one(1.0);
362 b->AddDomainIntegrator(new DomainLFIntegrator(one));
363 b->Assemble();
364
365 // 11. Define the solution vector x as a parallel finite element grid
366 // function corresponding to fespace. Initialize x with initial guess of
367 // zero, which satisfies the boundary conditions.
368 ParGridFunction x(fespace);
369 x = 0.0;
370
371 // 12. Set up the parallel bilinear form a(.,.) on the finite element space
372 // that will hold the matrix corresponding to the Laplacian operator.
373 ParBilinearForm *a = new ParBilinearForm(fespace);
374 ParBilinearForm *a_pc = NULL;
375 if (pc_choice == PCType::LOR) { a_pc = new ParBilinearForm(fespace_lor); }
376 if (pc_choice == PCType::HO) { a_pc = new ParBilinearForm(fespace); }
377
378 // 13. Assemble the parallel bilinear form and the corresponding linear
379 // system, applying any necessary transformations such as: parallel
380 // assembly, eliminating boundary conditions, applying conforming
381 // constraints for non-conforming AMR, static condensation, etc.
382 if (static_cond)
383 {
384 a->EnableStaticCondensation();
385 MFEM_VERIFY(pc_choice != PCType::LOR,
386 "cannot use LOR preconditioner with static condensation");
387 }
388
389 if (myid == 0)
390 {
391 cout << "Assembling the matrix ..." << flush;
392 }
393 tic_toc.Clear();
394 tic_toc.Start();
395 // Pre-allocate sparsity assuming dense element matrices
396 a->UsePrecomputedSparsity();
397
398 HPCBilinearForm *a_hpc = NULL;
399 Operator *a_oper = NULL;
400
401 if (!perf)
402 {
403 // Standard assembly using a diffusion domain integrator
404 a->AddDomainIntegrator(new DiffusionIntegrator(one));
405 a->Assemble();
406 }
407 else
408 {
409 // High-performance assembly/evaluation using the templated operator type
410 a_hpc = new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
411 if (matrix_free)
412 {
413 a_hpc->Assemble(); // partial assembly
414 }
415 else
416 {
417 a_hpc->AssembleBilinearForm(*a); // full matrix assembly
418 }
419 }
420 tic_toc.Stop();
421 if (myid == 0)
422 {
423 cout << " done, " << tic_toc.RealTime() << "s." << endl;
424 }
425
426 // 14. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
427 // preconditioner from hypre.
428
429 // Setup the operator matrix (if applicable)
431 Vector B, X;
432 if (perf && matrix_free)
433 {
434 a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
435 HYPRE_BigInt glob_size = fespace->GlobalTrueVSize();
436 if (myid == 0)
437 {
438 cout << "Size of linear system: " << glob_size << endl;
439 }
440 }
441 else
442 {
443 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
444 HYPRE_BigInt glob_size = A.GetGlobalNumRows();
445 if (myid == 0)
446 {
447 cout << "Size of linear system: " << glob_size << endl;
448 }
449 a_oper = &A;
450 }
451
452 // Setup the matrix used for preconditioning
453 if (myid == 0)
454 {
455 cout << "Assembling the preconditioning matrix ..." << flush;
456 }
457 tic_toc.Clear();
458 tic_toc.Start();
459
460 HypreParMatrix A_pc;
461 if (pc_choice == PCType::LOR)
462 {
463 // TODO: assemble the LOR matrix using the performance code
466 a_pc->Assemble();
467 a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
468 }
469 else if (pc_choice == PCType::HO)
470 {
471 if (!matrix_free)
472 {
473 A_pc.MakeRef(A); // matrix already assembled, reuse it
474 }
475 else
476 {
478 a_hpc->AssembleBilinearForm(*a_pc);
479 a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
480 }
481 }
482 tic_toc.Stop();
483 if (myid == 0)
484 {
485 cout << " done, " << tic_toc.RealTime() << "s." << endl;
486 }
487
488 // Solve with CG or PCG, depending if the matrix A_pc is available
489 CGSolver *pcg;
490 pcg = new CGSolver(MPI_COMM_WORLD);
491 pcg->SetRelTol(1e-6);
492 pcg->SetMaxIter(500);
493 pcg->SetPrintLevel(1);
494
495 HypreSolver *amg = NULL;
496
497 pcg->SetOperator(*a_oper);
498 if (pc_choice != PCType::NONE)
499 {
500 amg = new HypreBoomerAMG(A_pc);
501 pcg->SetPreconditioner(*amg);
502 }
503
504 tic_toc.Clear();
505 tic_toc.Start();
506
507 pcg->Mult(B, X);
508
509 tic_toc.Stop();
510 delete amg;
511
512 if (myid == 0)
513 {
514 // Note: In the pcg algorithm, the number of operator Mult() calls is
515 // N_iter and the number of preconditioner Mult() calls is N_iter+1.
516 cout << "Time per CG step: "
517 << tic_toc.RealTime() / pcg->GetNumIterations() << "s." << endl;
518 }
519
520 // 15. Recover the parallel grid function corresponding to X. This is the
521 // local finite element solution on each processor.
522 if (perf && matrix_free)
523 {
524 a_hpc->RecoverFEMSolution(X, *b, x);
525 }
526 else
527 {
528 a->RecoverFEMSolution(X, *b, x);
529 }
530
531 // 16. Save the refined mesh and the solution in parallel. This output can
532 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
533 {
534 ostringstream mesh_name, sol_name;
535 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
536 sol_name << "sol." << setfill('0') << setw(6) << myid;
537
538 ofstream mesh_ofs(mesh_name.str().c_str());
539 mesh_ofs.precision(8);
540 pmesh->Print(mesh_ofs);
541
542 ofstream sol_ofs(sol_name.str().c_str());
543 sol_ofs.precision(8);
544 x.Save(sol_ofs);
545 }
546
547 // 17. Send the solution by socket to a GLVis server.
548 if (visualization)
549 {
550 char vishost[] = "localhost";
551 socketstream sol_sock(vishost, visport);
552 sol_sock << "parallel " << num_procs << " " << myid << "\n";
553 sol_sock.precision(8);
554 sol_sock << "solution\n" << *pmesh << x << flush;
555 }
556
557 // 18. Free the used memory.
558 delete a;
559 delete a_hpc;
560 if (a_oper != &A) { delete a_oper; }
561 delete a_pc;
562 delete b;
563 delete fespace;
564 delete fespace_lor;
565 delete fec_lor;
566 if (order > 0) { delete fec; }
567 delete pmesh;
568 delete pcg;
569
570 return 0;
571}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:147
@ Positive
Bernstein polynomials.
Definition fe_base.hpp:37
@ ClosedUniform
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition fe_base.hpp:39
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition fe_base.hpp:96
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
Definition fe_base.hpp:115
void UsePrecomputedSparsity(int ps=1)
For scalar FE spaces, precompute the sparsity pattern of the matrix (assuming dense element matrices)...
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Conjugate gradient method.
Definition solvers.hpp:538
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:751
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
Definition solvers.hpp:551
A coefficient that is constant across space and time.
Class for domain integration .
Definition lininteg.hpp:106
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
virtual const char * Name() const
Definition fe_coll.hpp:79
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:275
The BoomerAMG solver in hypre.
Definition hypre.hpp:1717
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:408
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition hypre.hpp:699
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
Definition hypre.cpp:1474
Abstract class for hypre's solvers and preconditioners.
Definition hypre.hpp:1188
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.cpp:33
void SetRelTol(real_t rtol)
Definition solvers.hpp:229
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:174
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition solvers.hpp:280
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:72
void SetMaxIter(int max_it)
Definition solvers.hpp:231
Mesh data type.
Definition mesh.hpp:64
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:290
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:298
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1216
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9294
int MeshGenerator() const
Get the mesh generator/type.
Definition mesh.hpp:1243
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition mesh.cpp:6484
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
static int WorldRank()
Return the MPI rank in 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).
Abstract operator.
Definition operator.hpp:25
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.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A) override
Form the linear system matrix A, see FormLinearSystem() for details.
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:346
Class for parallel grid function.
Definition pgridfunc.hpp:50
Class for parallel linear form.
Class for parallel meshes.
Definition pmesh.hpp:34
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
Definition pmesh.cpp:1374
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4800
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition tic_toc.cpp:432
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition tic_toc.cpp:411
void Stop()
Stop the stopwatch.
Definition tic_toc.cpp:422
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
Definition tic_toc.cpp:406
Templated bilinear form class, cf. bilinearform.?pp.
The Integrator class combines a kernel and a coefficient.
Vector data type.
Definition vector.hpp:82
int dim
Definition ex24.cpp:53
int main()
HYPRE_Int HYPRE_BigInt
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
PCType
Definition ex1.cpp:48
const int mesh_p
Definition ex1p.cpp:57
const int sol_p
Definition ex1p.cpp:58
void mfem_error(const char *msg)
Definition error.cpp:154
StopWatch tic_toc
Definition tic_toc.cpp:450
const char vishost[]