MFEM  v4.1.0
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 Version
2 // PUMI Modification
3 //
4 // Compile with: make ex1p
5 //
6 // Sample runs:
7 // mpirun -np 8 ex1p -m ../../data/pumi/parallel/Kova/Kova100k_8.smb
8 // -p ../../data/pumi/geom/Kova.dmg -o 1 -go 2
9 //
10 // Note: Example models + meshes for the PUMI examples can be downloaded
11 // from github.com/mfem/data/pumi. After downloading we recommend
12 // creating a symbolic link to the above directory in ../../data.
13 //
14 // Description: This example code demonstrates the use of MFEM to define a
15 // simple finite element discretization of the Laplace problem
16 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
17 // Specifically, we discretize using a FE space of the specified
18 // order, or if order < 1 using an isoparametric/isogeometric
19 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
20 // NURBS mesh, etc.)
21 //
22 // The example highlights the use of mesh refinement, finite
23 // element grid functions, as well as linear and bilinear forms
24 // corresponding to the left-hand side and right-hand side of the
25 // discrete linear system. We also cover the explicit elimination
26 // of essential boundary conditions, static condensation, and the
27 // optional connection to the GLVis tool for visualization.
28 //
29 // This PUMI modification demonstrates how PUMI's API can be used
30 // to load a parallel PUMI mesh classified on a geometric model
31 // and then generate the corresponding parallel MFEM mesh. The
32 // example also performs a "uniform" refinement, similar to the
33 // MFEM examples, for coarse meshes. However, the refinement is
34 // performed using the PUMI API. The inputs are a Parasolid
35 // model, "*.xmt_txt" and SCOREC parallel meshes "*.smb". The
36 // option "-o" is used for the Finite Element order and "-go" for
37 // the geometry order. Note that they can be used independently:
38 // "-o 8 -go 3" solves for 8th order FE on third order geometry.
39 
40 #include "mfem.hpp"
41 #include <fstream>
42 #include <iostream>
43 
44 #ifdef MFEM_USE_SIMMETRIX
45 #include <SimUtil.h>
46 #include <gmi_sim.h>
47 #endif
48 #include <apfMDS.h>
49 #include <gmi_null.h>
50 #include <PCU.h>
51 #include <apfConvert.h>
52 #include <gmi_mesh.h>
53 #include <crv.h>
54 
55 using namespace std;
56 using namespace mfem;
57 
58 int main(int argc, char *argv[])
59 {
60  // 1. Initialize MPI.
61  int num_procs, myid;
62  MPI_Init(&argc, &argv);
63  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
64  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
65 
66  // 2. Parse command-line options.
67  const char *mesh_file = "../../data/pumi/parallel/Kova/Kova100k_8.smb";
68 #ifdef MFEM_USE_SIMMETRIX
69  const char *model_file = "../../data/pumi/geom/Kova.x_t";
70 #else
71  const char *model_file = "../../data/pumi/geom/Kova.dmg";
72 #endif
73  int order = 1;
74  bool static_cond = false;
75  bool visualization = 1;
76  int geom_order = 1;
77 
78  OptionsParser args(argc, argv);
79  args.AddOption(&mesh_file, "-m", "--mesh",
80  "Mesh file to use.");
81  args.AddOption(&order, "-o", "--order",
82  "Finite element order (polynomial degree) or -1 for"
83  " isoparametric space.");
84  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
85  "--no-static-condensation", "Enable static condensation.");
86  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
87  "--no-visualization",
88  "Enable or disable GLVis visualization.");
89  args.AddOption(&model_file, "-p", "--parasolid",
90  "Parasolid model to use.");
91  args.AddOption(&geom_order, "-go", "--geometry_order",
92  "Geometric order of the model");
93  args.Parse();
94  if (!args.Good())
95  {
96  if (myid == 0)
97  {
98  args.PrintUsage(cout);
99  }
100  MPI_Finalize();
101  return 1;
102  }
103  if (myid == 0)
104  {
105  args.PrintOptions(cout);
106  }
107 
108  // 3. Read the SCOREC Mesh
109  PCU_Comm_Init();
110 #ifdef MFEM_USE_SIMMETRIX
111  Sim_readLicenseFile(0);
112  gmi_sim_start();
113  gmi_register_sim();
114 #endif
115  gmi_register_mesh();
116 
117  apf::Mesh2* pumi_mesh;
118  pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
119 
120  // 4. Increase the geometry order and refine the mesh if necessary. Parallel
121  // uniform refinement is performed if the total number of elements is less
122  // than 10,000.
123  int dim = pumi_mesh->getDimension();
124  int nEle = pumi_mesh->count(dim);
125  int ref_levels = (int)floor(log(10000./nEle)/log(2.)/dim);
126 
127  if (geom_order > 1)
128  {
129  crv::BezierCurver bc(pumi_mesh, geom_order, 2);
130  bc.run();
131  }
132 
133  // Perform Uniform refinement
134  if (ref_levels > 1)
135  {
136  ma::Input* uniInput = ma::configureUniformRefine(pumi_mesh, ref_levels);
137 
138  if (geom_order > 1)
139  {
140  crv::adapt(uniInput);
141  }
142  else
143  {
144  ma::adapt(uniInput);
145  }
146  }
147 
148  pumi_mesh->verify();
149 
150  // 5. Create the parallel MFEM mesh object from the parallel PUMI mesh.
151  // We can handle triangular and tetrahedral meshes. Note that the
152  // mesh resolution is performed on the PUMI mesh.
153  ParMesh *pmesh = new ParPumiMesh(MPI_COMM_WORLD, pumi_mesh);
154 
155  // 6. Define a parallel finite element space on the parallel mesh. Here we
156  // use continuous Lagrange finite elements of the specified order. If
157  // order < 1, we instead use an isoparametric/isogeometric space.
159  if (order > 0)
160  {
161  fec = new H1_FECollection(order, dim);
162  }
163  else if (pmesh->GetNodes())
164  {
165  fec = pmesh->GetNodes()->OwnFEC();
166  if (myid == 0)
167  {
168  cout << "Using isoparametric FEs: " << fec->Name() << endl;
169  }
170  }
171  else
172  {
173  fec = new H1_FECollection(order = 1, dim);
174  }
175  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh, fec);
176  HYPRE_Int size = fespace->GlobalTrueVSize();
177  if (myid == 0)
178  {
179  cout << "Number of finite element unknowns: " << size << endl;
180  }
181 
182  // 7. Determine the list of true (i.e. parallel conforming) essential
183  // boundary dofs. In this example, the boundary conditions are defined
184  // by marking all the boundary attributes from the mesh as essential
185  // (Dirichlet) and converting them to a list of true dofs.
186  Array<int> ess_tdof_list;
187  if (pmesh->bdr_attributes.Size())
188  {
189  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
190  ess_bdr = 1;
191  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
192  }
193 
194  // 8. Set up the parallel linear form b(.) which corresponds to the
195  // right-hand side of the FEM linear system, which in this case is
196  // (1,phi_i) where phi_i are the basis functions in fespace.
197  ParLinearForm *b = new ParLinearForm(fespace);
198  ConstantCoefficient one(1.0);
200  b->Assemble();
201 
202  // 9. Define the solution vector x as a parallel finite element grid function
203  // corresponding to fespace. Initialize x with initial guess of zero,
204  // which satisfies the boundary conditions.
205  ParGridFunction x(fespace);
206  x = 0.0;
207 
208  // 10. Set up the parallel bilinear form a(.,.) on the finite element space
209  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
210  // domain integrator.
211  ParBilinearForm *a = new ParBilinearForm(fespace);
213 
214  // 11. Assemble the parallel bilinear form and the corresponding linear
215  // system, applying any necessary transformations such as: parallel
216  // assembly, eliminating boundary conditions, applying conforming
217  // constraints for non-conforming AMR, static condensation, etc.
218  if (static_cond) { a->EnableStaticCondensation(); }
219  a->Assemble();
220 
221  HypreParMatrix A;
222  Vector B, X;
223  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
224 
225  if (myid == 0)
226  {
227  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
228  }
229 
230  // 12. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
231  // preconditioner from hypre.
232  HypreSolver *amg = new HypreBoomerAMG(A);
233  HyprePCG *pcg = new HyprePCG(A);
234  pcg->SetTol(1e-12);
235  pcg->SetMaxIter(200);
236  pcg->SetPrintLevel(2);
237  pcg->SetPreconditioner(*amg);
238  pcg->Mult(B, X);
239 
240  // 13. Recover the parallel grid function corresponding to X. This is the
241  // local finite element solution on each processor.
242  a->RecoverFEMSolution(X, *b, x);
243 
244  // 14. Save the refined mesh and the solution in parallel. This output can
245  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
246  {
247  ostringstream mesh_name, sol_name;
248  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
249  sol_name << "sol." << setfill('0') << setw(6) << myid;
250 
251  ofstream mesh_ofs(mesh_name.str().c_str());
252  mesh_ofs.precision(8);
253  pmesh->Print(mesh_ofs);
254 
255  ofstream sol_ofs(sol_name.str().c_str());
256  sol_ofs.precision(8);
257  x.Save(sol_ofs);
258  }
259 
260  // 15. Send the solution by socket to a GLVis server.
261  if (visualization)
262  {
263  char vishost[] = "localhost";
264  int visport = 19916;
265  socketstream sol_sock(vishost, visport);
266  sol_sock << "parallel " << num_procs << " " << myid << "\n";
267  sol_sock.precision(8);
268  sol_sock << "solution\n" << *pmesh << x << flush;
269  }
270 
271  // 16. Free the used memory.
272  delete pcg;
273  delete amg;
274  delete a;
275  delete b;
276  delete fespace;
277  if (order > 0) { delete fec; }
278  delete pmesh;
279 
280  pumi_mesh->destroyNative();
281  apf::destroyMesh(pumi_mesh);
282  PCU_Comm_Free();
283 
284 #ifdef MFEM_USE_SIMMETRIX
285  gmi_sim_stop();
286  Sim_unregisterAllKeys();
287 #endif
288 
289  MPI_Finalize();
290 
291  return 0;
292 }
void SetTol(double tol)
Definition: hypre.cpp:2305
int Size() const
Logical size of the array.
Definition: array.hpp:124
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:770
Subclass constant coefficient.
Definition: coefficient.hpp:67
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()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:62
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2320
Class for PUMI parallel meshes.
Definition: pumi.hpp:72
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:421
double b
Definition: lissajous.cpp:42
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4102
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2310
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
PCG solver in hypre.
Definition: hypre.hpp:743
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:76
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:53
int dim
Definition: ex24.cpp:43
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2325
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
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:684
Vector data type.
Definition: vector.hpp:48
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6203
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2348
Class for parallel meshes.
Definition: pmesh.hpp:32
void EnableStaticCondensation()
bool Good() const
Definition: optparser.hpp:125