MFEM  v4.1.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 #include "mfem.hpp"
37 #include <fstream>
38 #include <iostream>
39 
40 #ifdef MFEM_USE_SIMMETRIX
41 #include <SimUtil.h>
42 #include <gmi_sim.h>
43 #endif
44 #include <apfMDS.h>
45 #include <gmi_null.h>
46 #include <PCU.h>
47 #include <apfConvert.h>
48 #include <gmi_mesh.h>
49 #include <crv.h>
50 
51 using namespace std;
52 using namespace mfem;
53 
54 int main(int argc, char *argv[])
55 {
56  // 1. Initialize MPI (required by PUMI).
57  int num_procs, myid;
58  MPI_Init(&argc, &argv);
59  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
60  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
61 
62  // 2. Parse command-line options.
63  const char *mesh_file = "../../data/pumi/serial/Kova.smb";
64 #ifdef MFEM_USE_SIMMETRIX
65  const char *model_file = "../../data/pumi/geom/Kova.x_t";
66 #else
67  const char *model_file = "../../data/pumi/geom/Kova.dmg";
68 #endif
69  int order = 1;
70  bool static_cond = false;
71  bool visualization = 1;
72  int geom_order = 1;
73 
74  OptionsParser args(argc, argv);
75  args.AddOption(&mesh_file, "-m", "--mesh",
76  "Mesh file to use.");
77  args.AddOption(&order, "-o", "--order",
78  "Finite element order (polynomial degree) or -1 for"
79  " isoparametric space.");
80  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
81  "--no-static-condensation", "Enable static condensation.");
82  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
83  "--no-visualization",
84  "Enable or disable GLVis visualization.");
85  args.AddOption(&model_file, "-p", "--parasolid",
86  "Parasolid model to use.");
87  args.AddOption(&geom_order, "-go", "--geometry_order",
88  "Geometric order of the model");
89  args.Parse();
90  if (!args.Good())
91  {
92  if (myid == 0)
93  {
94  args.PrintUsage(cout);
95  }
96  MPI_Finalize();
97  return 1;
98  }
99  if (myid == 0)
100  {
101  args.PrintOptions(cout);
102  }
103 
104  // 3. Read the SCOREC Mesh.
105  PCU_Comm_Init();
106 #ifdef MFEM_USE_SIMMETRIX
107  Sim_readLicenseFile(0);
108  gmi_sim_start();
109  gmi_register_sim();
110 #endif
111  gmi_register_mesh();
112 
113  apf::Mesh2* pumi_mesh;
114  pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
115 
116  // 4. Increase the geometry order if necessary.
117  if (geom_order > 1)
118  {
119  crv::BezierCurver bc(pumi_mesh, geom_order, 2);
120  bc.run();
121  }
122 
123  pumi_mesh->verify();
124 
125  // 5. Create the MFEM mesh object from the PUMI mesh. We can handle
126  // triangular and tetrahedral meshes. Other inputs are the same as the
127  // MFEM default constructor.
128  Mesh *mesh = new PumiMesh(pumi_mesh, 1, 1);
129  int dim = mesh->Dimension();
130 
131  // 6. Refine the mesh to increase the resolution. In this example we do
132  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
133  // largest number that gives a final mesh with no more than 50,000
134  // elements.
135  {
136  int ref_levels =
137  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
138  for (int l = 0; l < ref_levels; l++)
139  {
140  mesh->UniformRefinement();
141  }
142  }
143 
144  // 7. Define a finite element space on the mesh. Here we use continuous
145  // Lagrange finite elements of the specified order. If order < 1, we
146  // instead use an isoparametric/isogeometric space.
148  if (order > 0)
149  {
150  fec = new H1_FECollection(order, dim);
151  }
152  else if (mesh->GetNodes())
153  {
154  fec = mesh->GetNodes()->OwnFEC();
155  cout << "Using isoparametric FEs: " << fec->Name() << endl;
156  }
157  else
158  {
159  fec = new H1_FECollection(order = 1, dim);
160  }
161  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
162  cout << "Number of finite element unknowns: "
163  << fespace->GetTrueVSize() << endl;
164 
165  // 8. Determine the list of true (i.e. conforming) essential boundary dofs.
166  // In this example, the boundary conditions are defined by marking all
167  // the boundary attributes from the mesh as essential (Dirichlet) and
168  // converting them to a list of true dofs.
169  Array<int> ess_tdof_list;
170  if (mesh->bdr_attributes.Size())
171  {
172  Array<int> ess_bdr(mesh->bdr_attributes.Max());
173  ess_bdr = 1;
174  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
175  }
176 
177  // 9. Set up the linear form b(.) which corresponds to the right-hand side of
178  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
179  // the basis functions in the finite element fespace.
180  LinearForm *b = new LinearForm(fespace);
181  ConstantCoefficient one(1.0);
183  b->Assemble();
184 
185  // 10. Define the solution vector x as a finite element grid function
186  // corresponding to fespace. Initialize x with initial guess of zero,
187  // which satisfies the boundary conditions.
188  GridFunction x(fespace);
189  x = 0.0;
190 
191  // 11. Set up the bilinear form a(.,.) on the finite element space
192  // corresponding to the Laplacian operator -Delta, by adding the
193  // Diffusion domain integrator.
194  BilinearForm *a = new BilinearForm(fespace);
196 
197  // 12. Assemble the bilinear form and the corresponding linear system,
198  // applying any necessary transformations such as: eliminating boundary
199  // conditions, applying conforming constraints for non-conforming AMR,
200  // static condensation, etc.
201  if (static_cond) { a->EnableStaticCondensation(); }
202  a->Assemble();
203 
204  SparseMatrix A;
205  Vector B, X;
206  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
207 
208  cout << "Size of linear system: " << A.Height() << endl;
209 
210 #ifndef MFEM_USE_SUITESPARSE
211  // 13. Define a simple symmetric Gauss-Seidel preconditioner and use it to
212  // solve the system A X = B with PCG.
213  GSSmoother M(A);
214  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
215 #else
216  // 13. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
217  UMFPackSolver umf_solver;
218  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
219  umf_solver.SetOperator(A);
220  umf_solver.Mult(B, X);
221 #endif
222 
223  // 14. Recover the solution as a finite element grid function.
224  a->RecoverFEMSolution(X, *b, x);
225 
226  // 15. Save the refined mesh and the solution. This output can be viewed later
227  // using GLVis: "glvis -m refined.mesh -g sol.gf".
228  ofstream mesh_ofs("refined.mesh");
229  mesh_ofs.precision(8);
230  mesh->Print(mesh_ofs);
231  ofstream sol_ofs("sol.gf");
232  sol_ofs.precision(8);
233  x.Save(sol_ofs);
234 
235  // 16. Send the solution by socket to a GLVis server.
236  if (visualization)
237  {
238  char vishost[] = "localhost";
239  int visport = 19916;
240  socketstream sol_sock(vishost, visport);
241  sol_sock.precision(8);
242  sol_sock << "solution\n" << *mesh << x << flush;
243  }
244 
245  // 17. Free the used memory.
246  delete a;
247  delete b;
248  delete fespace;
249  if (order > 0) { delete fec; }
250  delete mesh;
251 
252  pumi_mesh->destroyNative();
253  apf::destroyMesh(pumi_mesh);
254  PCU_Comm_Free();
255 #ifdef MFEM_USE_SIMMETRIX
256  gmi_sim_stop();
257  Sim_unregisterAllKeys();
258 #endif
259 
260  MPI_Finalize();
261  return 0;
262 }
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 Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:67
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)
Definition: fespace.cpp:366
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
Base class for PUMI meshes.
Definition: pumi.hpp:44
int main(int argc, char *argv[])
Definition: ex1.cpp:62
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:580
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2731
Data type sparse matrix.
Definition: sparsemat.hpp:40
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
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:529
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
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:189
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:591
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
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
double a
Definition: lissajous.cpp:41
virtual const char * Name() const
Definition: fe_coll.hpp:53
int dim
Definition: ex24.cpp:43
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
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 linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:2382
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:2285
bool Good() const
Definition: optparser.hpp:125