MFEM v4.7.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 Laplace 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);
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 bool visualization = 1;
110
111 OptionsParser args(argc, argv);
112 args.AddOption(&mesh_file, "-m", "--mesh",
113 "Mesh file to use.");
114 args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
115 "Number of times to refine the mesh uniformly in serial;"
116 " -1 = auto: <= 10,000 elements.");
117 args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
118 "Number of times to refine the mesh uniformly in parallel.");
119 args.AddOption(&order, "-o", "--order",
120 "Finite element order (polynomial degree) or -1 for"
121 " isoparametric space.");
122 args.AddOption(&basis_type, "-b", "--basis-type",
123 "Basis: G - Gauss-Lobatto, P - Positive, U - Uniform");
124 args.AddOption(&perf, "-perf", "--hpc-version", "-std", "--standard-version",
125 "Enable high-performance, tensor-based, assembly/evaluation.");
126 args.AddOption(&matrix_free, "-mf", "--matrix-free", "-asm", "--assembly",
127 "Use matrix-free evaluation or efficient matrix assembly in "
128 "the high-performance version.");
129 args.AddOption(&pc, "-pc", "--preconditioner",
130 "Preconditioner: lor - low-order-refined (matrix-free) AMG, "
131 "ho - high-order (assembled) AMG, none.");
132 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
133 "--no-static-condensation", "Enable static condensation.");
134 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
135 "--no-visualization",
136 "Enable or disable GLVis visualization.");
137 args.Parse();
138 if (!args.Good())
139 {
140 if (myid == 0)
141 {
142 args.PrintUsage(cout);
143 }
144 return 1;
145 }
146 if (static_cond && perf && matrix_free)
147 {
148 if (myid == 0)
149 {
150 cout << "\nStatic condensation can not be used with matrix-free"
151 " evaluation!\n" << endl;
152 }
153 return 2;
154 }
155 MFEM_VERIFY(perf || !matrix_free,
156 "--standard-version is not compatible with --matrix-free");
157 if (myid == 0)
158 {
159 args.PrintOptions(cout);
160 }
161
162 PCType pc_choice;
163 if (!strcmp(pc, "ho")) { pc_choice = PCType::HO; }
164 else if (!strcmp(pc, "lor")) { pc_choice = PCType::LOR; }
165 else if (!strcmp(pc, "none")) { pc_choice = PCType::NONE; }
166 else
167 {
168 mfem_error("Invalid Preconditioner specified");
169 return 3;
170 }
171
172 if (myid == 0)
173 {
174 cout << "\nMFEM SIMD width: " << MFEM_SIMD_BYTES/sizeof(double)
175 << " doubles\n" << endl;
176 }
177
178 // See class BasisType in fem/fe_coll.hpp for available basis types
179 int basis = BasisType::GetType(basis_type[0]);
180 if (myid == 0)
181 {
182 cout << "Using " << BasisType::Name(basis) << " basis ..." << endl;
183 }
184
185 // 3. Read the (serial) mesh from the given mesh file on all processors. We
186 // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
187 // and volume meshes with the same code.
188 Mesh *mesh = new Mesh(mesh_file, 1, 1);
189 int dim = mesh->Dimension();
190
191 if (dim == 2)
192 {
193 return ex1_t<2>::run(mesh, ser_ref_levels, par_ref_levels, order, basis,
194 static_cond, pc_choice, perf, matrix_free,
195 visualization);
196 }
197 else if (dim == 3)
198 {
199 return ex1_t<3>::run(mesh, ser_ref_levels, par_ref_levels, order,
200 basis, static_cond, pc_choice, perf, matrix_free,
201 visualization);
202 }
203 else
204 {
205 MFEM_ABORT("Dimension must be 2 or 3.")
206 }
207
208 return 0;
209}
210
211template <int dim>
212int ex1_t<dim>::run(Mesh *mesh, int ser_ref_levels, int par_ref_levels,
213 int order, int basis, bool static_cond, PCType pc_choice,
214 bool perf, bool matrix_free, bool visualization)
215{
216 int num_procs, myid;
217 MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
218 MPI_Comm_rank(MPI_COMM_WORLD, &myid);
219
220 // 4. Check if the optimized version matches the given mesh
221 if (perf)
222 {
223 if (myid == 0)
224 {
225 cout << "High-performance version using integration rule with "
226 << int_rule_t::qpts << " points ..." << endl;
227 }
228 if (!mesh_t::MatchesGeometry(*mesh))
229 {
230 if (myid == 0)
231 {
232 cout << "The given mesh does not match the optimized 'geom' parameter.\n"
233 << "Recompile with suitable 'geom' value." << endl;
234 }
235 delete mesh;
236 return 4;
237 }
238 else if (!mesh_t::MatchesNodes(*mesh))
239 {
240 if (myid == 0)
241 {
242 cout << "Switching the mesh curvature to match the "
243 << "optimized value (order " << mesh_p << ") ..." << endl;
244 }
245 mesh->SetCurvature(mesh_p, false, -1, Ordering::byNODES);
246 }
247 }
248
249 // 5. Refine the serial mesh on all processors to increase the resolution. In
250 // this example we do 'ref_levels' of uniform refinement. We choose
251 // 'ref_levels' to be the largest number that gives a final mesh with no
252 // more than 10,000 elements, or as specified on the command line with the
253 // option '--refine-serial'.
254 {
255 int ref_levels = (ser_ref_levels != -1) ? ser_ref_levels :
256 (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
257 for (int l = 0; l < ref_levels; l++)
258 {
259 mesh->UniformRefinement();
260 }
261 }
262 if (!perf && mesh->NURBSext)
263 {
264 const int new_mesh_p = std::min(sol_p, mesh_p);
265 if (myid == 0)
266 {
267 cout << "NURBS mesh: switching the mesh curvature to be "
268 << "min(sol_p, mesh_p) = " << new_mesh_p << " ..." << endl;
269 }
270 mesh->SetCurvature(new_mesh_p, false, -1, Ordering::byNODES);
271 }
272
273 // 6. Define a parallel mesh by a partitioning of the serial mesh. Refine
274 // this mesh further in parallel to increase the resolution. Once the
275 // parallel mesh is defined, the serial mesh can be deleted.
276 ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
277 delete mesh;
278 {
279 for (int l = 0; l < par_ref_levels; l++)
280 {
281 pmesh->UniformRefinement();
282 }
283 }
284 if (pmesh->MeshGenerator() & 1) // simplex mesh
285 {
286 MFEM_VERIFY(pc_choice != PCType::LOR, "triangle and tet meshes do not "
287 "support the LOR preconditioner yet");
288 }
289
290 // 7. Define a parallel finite element space on the parallel mesh. Here we
291 // use continuous Lagrange finite elements of the specified order. If
292 // order < 1, we instead use an isoparametric/isogeometric space.
294 if (order > 0)
295 {
296 fec = new H1_FECollection(order, dim, basis);
297 }
298 else if (pmesh->GetNodes())
299 {
300 fec = pmesh->GetNodes()->OwnFEC();
301 if (myid == 0)
302 {
303 cout << "Using isoparametric FEs: " << fec->Name() << endl;
304 }
305 }
306 else
307 {
308 fec = new H1_FECollection(order = 1, dim, basis);
309 }
310 ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
311 HYPRE_BigInt size = fespace->GlobalTrueVSize();
312 if (myid == 0)
313 {
314 cout << "Number of finite element unknowns: " << size << endl;
315 }
316
317 ParMesh pmesh_lor;
318 FiniteElementCollection *fec_lor = NULL;
319 ParFiniteElementSpace *fespace_lor = NULL;
320 if (pc_choice == PCType::LOR)
321 {
322 int basis_lor = basis;
323 if (basis == BasisType::Positive) { basis_lor=BasisType::ClosedUniform; }
324 pmesh_lor = ParMesh::MakeRefined(*pmesh, order, basis_lor);
325 fec_lor = new H1_FECollection(1, dim);
326 fespace_lor = new ParFiniteElementSpace(&pmesh_lor, fec_lor);
327 }
328
329 // 8. Check if the optimized version matches the given space
330 if (perf && !sol_fes_t::Matches(*fespace))
331 {
332 if (myid == 0)
333 {
334 cout << "The given order does not match the optimized parameter.\n"
335 << "Recompile with suitable 'sol_p' value." << endl;
336 }
337 delete fespace;
338 delete fec;
339 delete mesh;
340 return 5;
341 }
342
343 // 9. Determine the list of true (i.e. parallel conforming) essential
344 // boundary dofs. In this example, the boundary conditions are defined
345 // by marking all the boundary attributes from the mesh as essential
346 // (Dirichlet) and converting them to a list of true dofs.
347 Array<int> ess_tdof_list;
348 if (pmesh->bdr_attributes.Size())
349 {
350 Array<int> ess_bdr(pmesh->bdr_attributes.Max());
351 ess_bdr = 1;
352 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
353 }
354
355 // 10. Set up the parallel linear form b(.) which corresponds to the
356 // right-hand side of the FEM linear system, which in this case is
357 // (1,phi_i) where phi_i are the basis functions in fespace.
358 ParLinearForm *b = new ParLinearForm(fespace);
359 ConstantCoefficient one(1.0);
360 b->AddDomainIntegrator(new DomainLFIntegrator(one));
361 b->Assemble();
362
363 // 11. Define the solution vector x as a parallel finite element grid
364 // function corresponding to fespace. Initialize x with initial guess of
365 // zero, which satisfies the boundary conditions.
366 ParGridFunction x(fespace);
367 x = 0.0;
368
369 // 12. Set up the parallel bilinear form a(.,.) on the finite element space
370 // that will hold the matrix corresponding to the Laplacian operator.
371 ParBilinearForm *a = new ParBilinearForm(fespace);
372 ParBilinearForm *a_pc = NULL;
373 if (pc_choice == PCType::LOR) { a_pc = new ParBilinearForm(fespace_lor); }
374 if (pc_choice == PCType::HO) { a_pc = new ParBilinearForm(fespace); }
375
376 // 13. Assemble the parallel bilinear form and the corresponding linear
377 // system, applying any necessary transformations such as: parallel
378 // assembly, eliminating boundary conditions, applying conforming
379 // constraints for non-conforming AMR, static condensation, etc.
380 if (static_cond)
381 {
382 a->EnableStaticCondensation();
383 MFEM_VERIFY(pc_choice != PCType::LOR,
384 "cannot use LOR preconditioner with static condensation");
385 }
386
387 if (myid == 0)
388 {
389 cout << "Assembling the matrix ..." << flush;
390 }
391 tic_toc.Clear();
392 tic_toc.Start();
393 // Pre-allocate sparsity assuming dense element matrices
394 a->UsePrecomputedSparsity();
395
396 HPCBilinearForm *a_hpc = NULL;
397 Operator *a_oper = NULL;
398
399 if (!perf)
400 {
401 // Standard assembly using a diffusion domain integrator
402 a->AddDomainIntegrator(new DiffusionIntegrator(one));
403 a->Assemble();
404 }
405 else
406 {
407 // High-performance assembly/evaluation using the templated operator type
408 a_hpc = new HPCBilinearForm(integ_t(coeff_t(1.0)), *fespace);
409 if (matrix_free)
410 {
411 a_hpc->Assemble(); // partial assembly
412 }
413 else
414 {
415 a_hpc->AssembleBilinearForm(*a); // full matrix assembly
416 }
417 }
418 tic_toc.Stop();
419 if (myid == 0)
420 {
421 cout << " done, " << tic_toc.RealTime() << "s." << endl;
422 }
423
424 // 14. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
425 // preconditioner from hypre.
426
427 // Setup the operator matrix (if applicable)
429 Vector B, X;
430 if (perf && matrix_free)
431 {
432 a_hpc->FormLinearSystem(ess_tdof_list, x, *b, a_oper, X, B);
433 HYPRE_BigInt glob_size = fespace->GlobalTrueVSize();
434 if (myid == 0)
435 {
436 cout << "Size of linear system: " << glob_size << endl;
437 }
438 }
439 else
440 {
441 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
442 HYPRE_BigInt glob_size = A.GetGlobalNumRows();
443 if (myid == 0)
444 {
445 cout << "Size of linear system: " << glob_size << endl;
446 }
447 a_oper = &A;
448 }
449
450 // Setup the matrix used for preconditioning
451 if (myid == 0)
452 {
453 cout << "Assembling the preconditioning matrix ..." << flush;
454 }
455 tic_toc.Clear();
456 tic_toc.Start();
457
458 HypreParMatrix A_pc;
459 if (pc_choice == PCType::LOR)
460 {
461 // TODO: assemble the LOR matrix using the performance code
464 a_pc->Assemble();
465 a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
466 }
467 else if (pc_choice == PCType::HO)
468 {
469 if (!matrix_free)
470 {
471 A_pc.MakeRef(A); // matrix already assembled, reuse it
472 }
473 else
474 {
476 a_hpc->AssembleBilinearForm(*a_pc);
477 a_pc->FormSystemMatrix(ess_tdof_list, A_pc);
478 }
479 }
480 tic_toc.Stop();
481 if (myid == 0)
482 {
483 cout << " done, " << tic_toc.RealTime() << "s." << endl;
484 }
485
486 // Solve with CG or PCG, depending if the matrix A_pc is available
487 CGSolver *pcg;
488 pcg = new CGSolver(MPI_COMM_WORLD);
489 pcg->SetRelTol(1e-6);
490 pcg->SetMaxIter(500);
491 pcg->SetPrintLevel(1);
492
493 HypreSolver *amg = NULL;
494
495 pcg->SetOperator(*a_oper);
496 if (pc_choice != PCType::NONE)
497 {
498 amg = new HypreBoomerAMG(A_pc);
499 pcg->SetPreconditioner(*amg);
500 }
501
502 tic_toc.Clear();
503 tic_toc.Start();
504
505 pcg->Mult(B, X);
506
507 tic_toc.Stop();
508 delete amg;
509
510 if (myid == 0)
511 {
512 // Note: In the pcg algorithm, the number of operator Mult() calls is
513 // N_iter and the number of preconditioner Mult() calls is N_iter+1.
514 cout << "Time per CG step: "
515 << tic_toc.RealTime() / pcg->GetNumIterations() << "s." << endl;
516 }
517
518 // 15. Recover the parallel grid function corresponding to X. This is the
519 // local finite element solution on each processor.
520 if (perf && matrix_free)
521 {
522 a_hpc->RecoverFEMSolution(X, *b, x);
523 }
524 else
525 {
526 a->RecoverFEMSolution(X, *b, x);
527 }
528
529 // 16. Save the refined mesh and the solution in parallel. This output can
530 // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
531 {
532 ostringstream mesh_name, sol_name;
533 mesh_name << "mesh." << setfill('0') << setw(6) << myid;
534 sol_name << "sol." << setfill('0') << setw(6) << myid;
535
536 ofstream mesh_ofs(mesh_name.str().c_str());
537 mesh_ofs.precision(8);
538 pmesh->Print(mesh_ofs);
539
540 ofstream sol_ofs(sol_name.str().c_str());
541 sol_ofs.precision(8);
542 x.Save(sol_ofs);
543 }
544
545 // 17. Send the solution by socket to a GLVis server.
546 if (visualization)
547 {
548 char vishost[] = "localhost";
549 int visport = 19916;
550 socketstream sol_sock(vishost, visport);
551 sol_sock << "parallel " << num_procs << " " << myid << "\n";
552 sol_sock.precision(8);
553 sol_sock << "solution\n" << *pmesh << x << flush;
554 }
555
556 // 18. Free the used memory.
557 delete a;
558 delete a_hpc;
559 if (a_oper != &A) { delete a_oper; }
560 delete a_pc;
561 delete b;
562 delete fespace;
563 delete fespace_lor;
564 delete fec_lor;
565 if (order > 0) { delete fec; }
566 delete pmesh;
567 delete pcg;
568
569 return 0;
570}
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:144
@ Positive
Bernstein polynomials.
Definition fe_base.hpp:33
@ ClosedUniform
Nodes: x_i = i/(n-1), i=0,...,n-1.
Definition fe_base.hpp:35
static const char * Name(int b_type)
Check and convert a BasisType constant to a string identifier.
Definition fe_base.hpp:92
static int GetType(char b_ident)
Convert char basis identifier to a BasisType constant.
Definition fe_base.hpp:111
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: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.
Class for domain integration .
Definition lininteg.hpp:109
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:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition hypre.hpp:679
void MakeRef(const HypreParMatrix &master)
Make this HypreParMatrix a reference to 'master'.
Definition hypre.cpp:1408
Abstract class for hypre's solvers and preconditioners.
Definition hypre.hpp:1162
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
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition solvers.hpp:260
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
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
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
int MeshGenerator() const
Get the mesh generator/type.
Definition mesh.hpp:1187
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:6211
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
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.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
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:285
Class for parallel grid function.
Definition pgridfunc.hpp:33
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:1375
void Print(std::ostream &out=mfem::out, const std::string &comments="") const override
Definition pmesh.cpp:4801
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:80
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
const int visport
void mfem_error(const char *msg)
Definition error.cpp:154
StopWatch tic_toc
Definition tic_toc.cpp:450
const char vishost[]