MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex1p.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - Parallel NURBS Version
2 //
3 // Compile with: make ex1p
4 //
5 // Sample runs: mpirun -np 4 ex1p -m ../../data/square-disc.mesh
6 // mpirun -np 4 ex1p -m ../../data/star.mesh
7 // mpirun -np 4 ex1p -m ../../data/escher.mesh
8 // mpirun -np 4 ex1p -m ../../data/fichera.mesh
9 // mpirun -np 4 ex1p -m ../../data/square-disc-p2.vtk -o 2
10 // mpirun -np 4 ex1p -m ../../data/square-disc-p3.mesh -o 3
11 // mpirun -np 4 ex1p -m ../../data/square-disc-nurbs.mesh -o -1
12 // mpirun -np 4 ex1p -m ../../data/disc-nurbs.mesh -o -1
13 // mpirun -np 4 ex1p -m ../../data/pipe-nurbs.mesh -o -1
14 // mpirun -np 4 ex1p -m ../../data/ball-nurbs.mesh -o 2
15 // mpirun -np 4 ex1p -m ../../data/star-surf.mesh
16 // mpirun -np 4 ex1p -m ../../data/square-disc-surf.mesh
17 // mpirun -np 4 ex1p -m ../../data/inline-segment.mesh
18 // mpirun -np 4 ex1p -m ../../data/amr-quad.mesh
19 // mpirun -np 4 ex1p -m ../../data/amr-hex.mesh
20 // mpirun -np 4 ex1p -m ../../data/mobius-strip.mesh
21 // mpirun -np 4 ex1p -m ../../data/mobius-strip.mesh -o -1 -sc
22 //
23 // Description: This example code demonstrates the use of MFEM to define a
24 // simple finite element discretization of the Laplace problem
25 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
26 // Specifically, we discretize using a FE space of the specified
27 // order, or if order < 1 using an isoparametric/isogeometric
28 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
29 // NURBS mesh, etc.)
30 //
31 // The example highlights the use of mesh refinement, finite
32 // element grid functions, as well as linear and bilinear forms
33 // corresponding to the left-hand side and right-hand side of the
34 // discrete linear system. We also cover the explicit elimination
35 // of essential boundary conditions, static condensation, and the
36 // optional connection to the GLVis tool for visualization.
37 
38 #include "mfem.hpp"
39 #include <fstream>
40 #include <iostream>
41 
42 using namespace std;
43 using namespace mfem;
44 
45 int main(int argc, char *argv[])
46 {
47  // 1. Initialize MPI.
48  int num_procs, myid;
49  MPI_Init(&argc, &argv);
50  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
51  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
52 
53  // 2. Parse command-line options.
54  const char *mesh_file = "../../data/star.mesh";
55  Array<int> order(1);
56  order[0] = 1;
57  bool static_cond = false;
58  bool visualization = 1;
59 
60  OptionsParser args(argc, argv);
61  args.AddOption(&mesh_file, "-m", "--mesh",
62  "Mesh file to use.");
63  args.AddOption(&order, "-o", "--order",
64  "Finite element order (polynomial degree) or -1 for"
65  " isoparametric space.");
66  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
67  "--no-static-condensation", "Enable static condensation.");
68  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
69  "--no-visualization",
70  "Enable or disable GLVis visualization.");
71  args.Parse();
72  if (!args.Good())
73  {
74  if (myid == 0)
75  {
76  args.PrintUsage(cout);
77  }
78  MPI_Finalize();
79  return 1;
80  }
81  if (myid == 0)
82  {
83  args.PrintOptions(cout);
84  }
85 
86  // 3. Read the (serial) mesh from the given mesh file on all processors. We
87  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
88  // and volume meshes with the same code.
89  Mesh *mesh = new Mesh(mesh_file, 1, 1);
90  int dim = mesh->Dimension();
91 
92  // 4. Refine the serial mesh on all processors to increase the resolution. In
93  // this example we do 'ref_levels' of uniform refinement. We choose
94  // 'ref_levels' to be the largest number that gives a final mesh with no
95  // more than 10,000 elements.
96  {
97  int ref_levels =
98  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
99  for (int l = 0; l < ref_levels; l++)
100  {
101  mesh->UniformRefinement();
102  }
103  }
104 
105  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
106  // this mesh further in parallel to increase the resolution. Once the
107  // parallel mesh is defined, the serial mesh can be deleted.
108  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
109  delete mesh;
110  {
111  int par_ref_levels = 2;
112  for (int l = 0; l < par_ref_levels; l++)
113  {
114  pmesh->UniformRefinement();
115  }
116  }
117 
118  // 6. Define a parallel finite element space on the parallel mesh. Here we
119  // use continuous Lagrange finite elements of the specified order. If
120  // order < 1, we instead use an isoparametric/isogeometric space.
122  NURBSExtension *NURBSext = NULL;
123  int own_fec = 0;
124 
125  if (order[0] == -1) // Isoparametric
126  {
127  if (pmesh->GetNodes())
128  {
129  fec = pmesh->GetNodes()->OwnFEC();
130  own_fec = 0;
131  cout << "Using isoparametric FEs: " << fec->Name() << endl;
132  }
133  else
134  {
135  cout <<"Mesh does not have FEs --> Assume order 1.\n";
136  fec = new H1_FECollection(1, dim);
137  own_fec = 1;
138  }
139  }
140  else if (pmesh->NURBSext && (order[0] > 0) ) // Subparametric NURBS
141  {
142  fec = new NURBSFECollection(order[0]);
143  own_fec = 1;
144  int nkv = pmesh->NURBSext->GetNKV();
145 
146  if (order.Size() == 1)
147  {
148  int tmp = order[0];
149  order.SetSize(nkv);
150  order = tmp;
151  }
152  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
153  NURBSext = new NURBSExtension(pmesh->NURBSext, order);
154  }
155  else
156  {
157  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
158  fec = new H1_FECollection(abs(order[0]), dim);
159  own_fec = 1;
160  }
161  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh,NURBSext,fec);
162  HYPRE_Int size = fespace->GlobalTrueVSize();
163  if (myid == 0)
164  {
165  cout << "Number of finite element unknowns: " << size << endl;
166  }
167 
168  // 7. Determine the list of true (i.e. parallel conforming) essential
169  // boundary dofs. In this example, the boundary conditions are defined
170  // by marking all the boundary attributes from the mesh as essential
171  // (Dirichlet) and converting them to a list of true dofs.
172  Array<int> ess_tdof_list;
173  if (pmesh->bdr_attributes.Size())
174  {
175  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
176  ess_bdr = 1;
177  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
178  }
179 
180  // 8. Set up the parallel linear form b(.) which corresponds to the
181  // right-hand side of the FEM linear system, which in this case is
182  // (1,phi_i) where phi_i are the basis functions in fespace.
183  ParLinearForm *b = new ParLinearForm(fespace);
184  ConstantCoefficient one(1.0);
186  b->Assemble();
187 
188  // 9. Define the solution vector x as a parallel finite element grid function
189  // corresponding to fespace. Initialize x with initial guess of zero,
190  // which satisfies the boundary conditions.
191  ParGridFunction x(fespace);
192  x = 0.0;
193 
194  // 10. Set up the parallel bilinear form a(.,.) on the finite element space
195  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
196  // domain integrator.
197  ParBilinearForm *a = new ParBilinearForm(fespace);
199 
200  // 11. Assemble the parallel bilinear form and the corresponding linear
201  // system, applying any necessary transformations such as: parallel
202  // assembly, eliminating boundary conditions, applying conforming
203  // constraints for non-conforming AMR, static condensation, etc.
204  if (static_cond) { a->EnableStaticCondensation(); }
205  a->Assemble();
206 
207  HypreParMatrix A;
208  Vector B, X;
209  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
210 
211  if (myid == 0)
212  {
213  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
214  }
215 
216  // 12. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
217  // preconditioner from hypre.
218  HypreSolver *amg = new HypreBoomerAMG(A);
219  HyprePCG *pcg = new HyprePCG(A);
220  pcg->SetTol(1e-12);
221  pcg->SetMaxIter(200);
222  pcg->SetPrintLevel(2);
223  pcg->SetPreconditioner(*amg);
224  pcg->Mult(B, X);
225 
226  // 13. Recover the parallel grid function corresponding to X. This is the
227  // local finite element solution on each processor.
228  a->RecoverFEMSolution(X, *b, x);
229 
230  // 14. Save the refined mesh and the solution in parallel. This output can
231  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
232  {
233  ostringstream mesh_name, sol_name;
234  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
235  sol_name << "sol." << setfill('0') << setw(6) << myid;
236 
237  ofstream mesh_ofs(mesh_name.str().c_str());
238  mesh_ofs.precision(8);
239  pmesh->Print(mesh_ofs);
240 
241  ofstream sol_ofs(sol_name.str().c_str());
242  sol_ofs.precision(8);
243  x.Save(sol_ofs);
244  }
245 
246  // 15. Send the solution by socket to a GLVis server.
247  if (visualization)
248  {
249  char vishost[] = "localhost";
250  int visport = 19916;
251  socketstream sol_sock(vishost, visport);
252  sol_sock << "parallel " << num_procs << " " << myid << "\n";
253  sol_sock.precision(8);
254  sol_sock << "solution\n" << *pmesh << x << flush;
255  }
256 
257  // 16. Save data in the VisIt format
258  VisItDataCollection visit_dc("Example1-Parallel", pmesh);
259  visit_dc.RegisterField("solution", &x);
260  visit_dc.Save();
261 
262  // 17. Free the used memory.
263  delete pcg;
264  delete amg;
265  delete a;
266  delete b;
267  delete fespace;
268  if (own_fec) { delete fec; }
269  delete pmesh;
270 
271  MPI_Finalize();
272 
273  return 0;
274 }
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:266
void SetTol(double tol)
Definition: hypre.cpp:2170
int Size() const
Logical size of the array.
Definition: array.hpp:133
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:626
Subclass constant coefficient.
Definition: coefficient.hpp:57
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
int GetNKV() const
Definition: nurbs.hpp:319
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:398
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:45
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2185
The BoomerAMG solver in hypre.
Definition: hypre.hpp:796
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Definition: hypre.hpp:379
int dim
Definition: ex3.cpp:47
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:146
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
Data collection with VisIt I/O routines.
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:247
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3428
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2175
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
PCG solver in hypre.
Definition: hypre.hpp:656
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)
Definition: optparser.hpp:74
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:569
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:176
virtual const char * Name() const
Definition: fe_coll.hpp:50
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2190
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:620
Vector data type.
Definition: vector.hpp:48
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5422
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:175
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2211
Class for parallel meshes.
Definition: pmesh.hpp:32
void EnableStaticCondensation()
bool Good() const
Definition: optparser.hpp:120