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