MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex1.cpp
Go to the documentation of this file.
1 // MFEM Example 1
2 // PUMI Modification
3 //
4 // Compile with: make ex1
5 //
6 // Sample runs:
7 // ex1 -m ../../data/pumi/serial/Kova.smb -p ../../data/pumi/geom/Kova.dmg
8 //
9 // Note: Example models + meshes for the PUMI examples can be downloaded
10 // from github.com/mfem/data/pumi. After downloading we recommend
11 // creating a symbolic link to the above directory in ../../data.
12 //
13 // Description: This example code demonstrates the use of MFEM to define a
14 // simple finite element discretization of the Laplace problem
15 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
16 // Specifically, we discretize using a FE space of the specified
17 // order, or if order < 1 using an isoparametric/isogeometric
18 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
19 // NURBS mesh, etc.)
20 //
21 // The example highlights the use of mesh refinement, finite
22 // element grid functions, as well as linear and bilinear forms
23 // corresponding to the left-hand side and right-hand side of the
24 // discrete linear system. We also cover the explicit elimination
25 // of essential boundary conditions, static condensation, and the
26 // optional connection to the GLVis tool for visualization.
27 //
28 // This PUMI modification demonstrates how PUMI's API can be used
29 // to load a PUMI mesh classified on a geometric model and then
30 // convert it to the MFEM mesh format. The inputs are a Parasolid
31 // model, "*.xmt_txt" and a SCOREC mesh "*.smb". The option "-o"
32 // is used for the Finite Element order and "-go" is used for the
33 // geometry order. Note that they can be used independently, i.e.
34 // "-o 8 -go 3" solves for 8th order FE on a third order geometry.
35 //
36 // NOTE: Model/Mesh files for this example are in the (large) data file
37 // repository of MFEM here https://github.com/mfem/data under the
38 // folder named "pumi", which consists of the following sub-folders:
39 // a) geom --> model files
40 // b) parallel --> parallel pumi mesh files
41 // c) serial --> serial pumi mesh files
42 
43 #include "mfem.hpp"
44 #include <fstream>
45 #include <iostream>
46 
47 #ifdef MFEM_USE_SIMMETRIX
48 #include <SimUtil.h>
49 #include <gmi_sim.h>
50 #endif
51 #include <apfMDS.h>
52 #include <gmi_null.h>
53 #include <PCU.h>
54 #include <apfConvert.h>
55 #include <gmi_mesh.h>
56 #include <crv.h>
57 
58 #ifndef MFEM_USE_PUMI
59 #error This example requires that MFEM is built with MFEM_USE_PUMI=YES
60 #endif
61 
62 using namespace std;
63 using namespace mfem;
64 
65 int main(int argc, char *argv[])
66 {
67  // 1. Initialize MPI (required by PUMI).
68  int num_procs, myid;
69  MPI_Init(&argc, &argv);
70  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
71  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
72 
73  // 2. Parse command-line options.
74  const char *mesh_file = "../../data/pumi/serial/Kova.smb";
75 #ifdef MFEM_USE_SIMMETRIX
76  const char *model_file = "../../data/pumi/geom/Kova.x_t";
77 #else
78  const char *model_file = "../../data/pumi/geom/Kova.dmg";
79 #endif
80  int order = 1;
81  bool static_cond = false;
82  bool visualization = 1;
83  int geom_order = 1;
84 
85  OptionsParser args(argc, argv);
86  args.AddOption(&mesh_file, "-m", "--mesh",
87  "Mesh file to use.");
88  args.AddOption(&order, "-o", "--order",
89  "Finite element order (polynomial degree) or -1 for"
90  " isoparametric space.");
91  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
92  "--no-static-condensation", "Enable static condensation.");
93  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
94  "--no-visualization",
95  "Enable or disable GLVis visualization.");
96  args.AddOption(&model_file, "-p", "--parasolid",
97  "Parasolid model to use.");
98  args.AddOption(&geom_order, "-go", "--geometry_order",
99  "Geometric order of the model");
100  args.Parse();
101  if (!args.Good())
102  {
103  if (myid == 0)
104  {
105  args.PrintUsage(cout);
106  }
107  MPI_Finalize();
108  return 1;
109  }
110  if (myid == 0)
111  {
112  args.PrintOptions(cout);
113  }
114 
115  // 3. Read the SCOREC Mesh.
116  PCU_Comm_Init();
117 #ifdef MFEM_USE_SIMMETRIX
118  Sim_readLicenseFile(0);
119  gmi_sim_start();
120  gmi_register_sim();
121 #endif
122  gmi_register_mesh();
123 
124  apf::Mesh2* pumi_mesh;
125  pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
126 
127  // 4. Increase the geometry order if necessary.
128  if (geom_order > 1)
129  {
130  crv::BezierCurver bc(pumi_mesh, geom_order, 2);
131  bc.run();
132  }
133 
134  pumi_mesh->verify();
135 
136  // 5. Create the MFEM mesh object from the PUMI mesh. We can handle
137  // triangular and tetrahedral meshes. Other inputs are the same as the
138  // MFEM default constructor.
139  Mesh *mesh = new PumiMesh(pumi_mesh, 1, 1);
140  int dim = mesh->Dimension();
141 
142  // 6. Refine the mesh to increase the resolution. In this example we do
143  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
144  // largest number that gives a final mesh with no more than 50,000
145  // elements.
146  {
147  int ref_levels =
148  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
149  for (int l = 0; l < ref_levels; l++)
150  {
151  mesh->UniformRefinement();
152  }
153  }
154 
155  // 7. Define a finite element space on the mesh. Here we use continuous
156  // Lagrange finite elements of the specified order. If order < 1, we
157  // instead use an isoparametric/isogeometric space.
159  if (order > 0)
160  {
161  fec = new H1_FECollection(order, dim);
162  }
163  else if (mesh->GetNodes())
164  {
165  fec = mesh->GetNodes()->OwnFEC();
166  cout << "Using isoparametric FEs: " << fec->Name() << endl;
167  }
168  else
169  {
170  fec = new H1_FECollection(order = 1, dim);
171  }
172  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
173  cout << "Number of finite element unknowns: "
174  << fespace->GetTrueVSize() << endl;
175 
176  // 8. Determine the list of true (i.e. conforming) essential boundary dofs.
177  // In this example, the boundary conditions are defined by marking all
178  // the boundary attributes from the mesh as essential (Dirichlet) and
179  // converting them to a list of true dofs.
180  Array<int> ess_tdof_list;
181  if (mesh->bdr_attributes.Size())
182  {
183  Array<int> ess_bdr(mesh->bdr_attributes.Max());
184  ess_bdr = 1;
185  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
186  }
187 
188  // 9. Set up the linear form b(.) which corresponds to the right-hand side of
189  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
190  // the basis functions in the finite element fespace.
191  LinearForm *b = new LinearForm(fespace);
192  ConstantCoefficient one(1.0);
194  b->Assemble();
195 
196  // 10. Define the solution vector x as a finite element grid function
197  // corresponding to fespace. Initialize x with initial guess of zero,
198  // which satisfies the boundary conditions.
199  GridFunction x(fespace);
200  x = 0.0;
201 
202  // 11. Set up the bilinear form a(.,.) on the finite element space
203  // corresponding to the Laplacian operator -Delta, by adding the
204  // Diffusion domain integrator.
205  BilinearForm *a = new BilinearForm(fespace);
207 
208  // 12. Assemble the bilinear form and the corresponding linear system,
209  // applying any necessary transformations such as: eliminating boundary
210  // conditions, applying conforming constraints for non-conforming AMR,
211  // static condensation, etc.
212  if (static_cond) { a->EnableStaticCondensation(); }
213  a->Assemble();
214 
215  SparseMatrix A;
216  Vector B, X;
217  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
218 
219  cout << "Size of linear system: " << A.Height() << endl;
220 
221 #ifndef MFEM_USE_SUITESPARSE
222  // 13. Define a simple symmetric Gauss-Seidel preconditioner and use it to
223  // solve the system A X = B with PCG.
224  GSSmoother M(A);
225  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
226 #else
227  // 13. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
228  UMFPackSolver umf_solver;
229  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
230  umf_solver.SetOperator(A);
231  umf_solver.Mult(B, X);
232 #endif
233 
234  // 14. Recover the solution as a finite element grid function.
235  a->RecoverFEMSolution(X, *b, x);
236 
237  // 15. Save the refined mesh and the solution. This output can be viewed later
238  // using GLVis: "glvis -m refined.mesh -g sol.gf".
239  ofstream mesh_ofs("refined.mesh");
240  mesh_ofs.precision(8);
241  mesh->Print(mesh_ofs);
242  ofstream sol_ofs("sol.gf");
243  sol_ofs.precision(8);
244  x.Save(sol_ofs);
245 
246  // 16. 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.precision(8);
253  sol_sock << "solution\n" << *mesh << x << flush;
254  }
255 
256  // 17. Free the used memory.
257  delete a;
258  delete b;
259  delete fespace;
260  if (order > 0) { delete fec; }
261  delete mesh;
262 
263  pumi_mesh->destroyNative();
264  apf::destroyMesh(pumi_mesh);
265  PCU_Comm_Free();
266 #ifdef MFEM_USE_SIMMETRIX
267  gmi_sim_stop();
268  Sim_unregisterAllKeys();
269 #endif
270 
271  MPI_Finalize();
272  return 0;
273 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1234
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
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(...
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:415
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:737
Base class for PUMI meshes.
Definition: pumi.hpp:45
int main(int argc, char *argv[])
Definition: ex1.cpp:66
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:731
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3417
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:65
constexpr char vishost[]
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:8382
constexpr int visport
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:733
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:403
int Dimension() const
Definition: mesh.hpp:788
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:434
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:201
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:742
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
double a
Definition: lissajous.cpp:41
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
virtual const char * Name() const
Definition: fe_coll.hpp:61
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:304
Vector data type.
Definition: vector.hpp:51
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:6603
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2818
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:2721
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:145