MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex2.cpp
Go to the documentation of this file.
1// MFEM Example 2
2// PUMI Modification
3//
4// Compile with: make ex2
5//
6// Sample runs:
7// ex2 -m ../../data/pumi/serial/pillbox.smb -p ../../data/pumi/geom/pillbox.dmg
8// -bf ../../data/pumi/serial/boundary.mesh
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 solves a simple linear elasticity problem
15// describing a multi-material cantilever beam.
16//
17// Specifically, we approximate the weak form of -div(sigma(u))=0
18// where sigma(u)=lambda*div(u)*I+mu*(grad*u+u*grad) is the stress
19// tensor corresponding to displacement field u, and lambda and mu
20// are the material Lame constants. The boundary conditions are
21// u=0 on the fixed part of the boundary with attribute 1, and
22// sigma(u).n=f on the remainder with f being a constant pull down
23// vector on boundary elements with attribute 2, and zero
24// otherwise. The geometry of the domain is assumed to be as
25// follows:
26// boundary
27// attribute 2
28// (push down)
29// ||
30// \/
31// +----------+
32// | |
33// | |
34// +---------| material |----------+
35// boundary --->| material| 2 | material |<--- boundary
36// attribute 1 | 1 | | 3 | attribute 1
37// (fixed) +---------+----------+----------+ (fixed)
38//
39// The example demonstrates the use of high-order and NURBS vector
40// finite element spaces with the linear elasticity bilinear form,
41// meshes with curved elements, and the definition of piece-wise
42// constant and vector coefficient objects. Static condensation is
43// also illustrated.
44//
45// We recommend viewing Example 1 before viewing this example.
46//
47// NOTE: Model/Mesh files for this example are in the (large) data file
48// repository of MFEM here https://github.com/mfem/data under the
49// folder named "pumi", which consists of the following sub-folders:
50// a) geom --> model files
51// b) parallel --> parallel pumi mesh files
52// c) serial --> serial pumi mesh files
53
54
55#include "mfem.hpp"
56#include <fstream>
57#include <iostream>
58
60
61#ifdef MFEM_USE_SIMMETRIX
62#include <SimUtil.h>
63#include <gmi_sim.h>
64#endif
65#include <apfMDS.h>
66#include <gmi_null.h>
67#include <PCU.h>
68#include <apfConvert.h>
69#include <gmi_mesh.h>
70#include <crv.h>
71
72#ifndef MFEM_USE_PUMI
73#error This example requires that MFEM is built with MFEM_USE_PUMI=YES
74#endif
75
76using namespace std;
77using namespace mfem;
78
79int main(int argc, char *argv[])
80{
81 // 1. Initialize MPI (required by PUMI) and HYPRE.
82 Mpi::Init(argc, argv);
84
85 // 2. Parse command-line options.
86 const char *mesh_file = "../../data/pumi/serial/pillbox.smb";
87 const char *boundary_file = "../../data/pumi/serial/boundary.mesh";
88#ifdef MFEM_USE_SIMMETRIX
89 const char *model_file = "../../data/pumi/geom/pillbox.smd";
90#else
91 const char *model_file = "../../data/pumi/geom/pillbox.dmg";
92#endif
93 int order = 1;
94 bool static_cond = false;
95 bool visualization = 1;
96 int geom_order = 1;
97
98 OptionsParser args(argc, argv);
99 args.AddOption(&mesh_file, "-m", "--mesh",
100 "Mesh file to use.");
101 args.AddOption(&order, "-o", "--order",
102 "Finite element order (polynomial degree).");
103 args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
104 "--no-static-condensation", "Enable static condensation.");
105 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
106 "--no-visualization",
107 "Enable or disable GLVis visualization.");
108 args.AddOption(&model_file, "-p", "--parasolid",
109 "Parasolid model to use.");
110 args.AddOption(&geom_order, "-go", "--geometry_order",
111 "Geometric order of the model");
112 args.AddOption(&boundary_file, "-bf", "--txt",
113 "txt file containing boundary tags");
114 args.Parse();
115 if (!args.Good())
116 {
117 args.PrintUsage(cout);
118 return 1;
119 }
120 args.PrintOptions(cout);
121
122 // 3. Read the SCOREC Mesh.
123 PCU_Comm_Init();
124#ifdef MFEM_USE_SIMMETRIX
125 Sim_readLicenseFile(0);
126 gmi_sim_start();
127 gmi_register_sim();
128#endif
129 gmi_register_mesh();
130
131 apf::Mesh2* pumi_mesh;
132 pumi_mesh = apf::loadMdsMesh(model_file, mesh_file);
133
134 // 4. Increase the geometry order if necessary.
135 if (geom_order > 1)
136 {
137 crv::BezierCurver bc(pumi_mesh, geom_order, 0);
138 bc.run();
139 }
140 pumi_mesh->verify();
141
142 // Read boundary
143 string bdr_tags;
144 named_ifgzstream input_bdr(boundary_file);
145 input_bdr >> ws;
146 getline(input_bdr, bdr_tags);
147 filter_dos(bdr_tags);
148 cout << " the boundary tag is : " << bdr_tags << endl;
149 Array<int> Dirichlet;
150 int numOfent;
151 if (bdr_tags == "Dirichlet")
152 {
153 input_bdr >> numOfent;
154 cout << " num of Dirichlet bdr conditions : " << numOfent << endl;
155 Dirichlet.SetSize(numOfent);
156 for (int kk = 0; kk < numOfent; kk++)
157 {
158 input_bdr >> Dirichlet[kk];
159 }
160 }
161 Dirichlet.Print();
162
163 Array<int> load_bdr;
164 skip_comment_lines(input_bdr, '#');
165 input_bdr >> bdr_tags;
166 filter_dos(bdr_tags);
167 cout << " the boundary tag is : " << bdr_tags << endl;
168 if (bdr_tags == "Load")
169 {
170 input_bdr >> numOfent;
171 load_bdr.SetSize(numOfent);
172 cout << " num of load bdr conditions : " << numOfent << endl;
173 for (int kk = 0; kk < numOfent; kk++)
174 {
175 input_bdr >> load_bdr[kk];
176 }
177 }
178 load_bdr.Print();
179
180 // 5. Create the MFEM mesh object from the PUMI mesh. We can handle triangular
181 // and tetrahedral meshes. Other inputs are the same as MFEM default
182 // constructor.
183 Mesh *mesh = new PumiMesh(pumi_mesh, 1, 1);
184 int dim = mesh->Dimension();
185
186 // Boundary conditions hack.
187 apf::MeshIterator* itr = pumi_mesh->begin(dim-1);
188 apf::MeshEntity* ent ;
189 int bdr_cnt = 0;
190 while ((ent = pumi_mesh->iterate(itr)))
191 {
192 apf::ModelEntity *me = pumi_mesh->toModel(ent);
193 if (pumi_mesh->getModelType(me) == (dim-1))
194 {
195 // Everywhere 3 as initial
196 (mesh->GetBdrElement(bdr_cnt))->SetAttribute(3);
197 int tag = pumi_mesh->getModelTag(me);
198 if (Dirichlet.Find(tag) != -1)
199 {
200 // Dirichlet attr -> 1
201 (mesh->GetBdrElement(bdr_cnt))->SetAttribute(1);
202 }
203 else if (load_bdr.Find(tag) != -1)
204 {
205 // Load attr -> 2
206 (mesh->GetBdrElement(bdr_cnt))->SetAttribute(2);
207 }
208 bdr_cnt++;
209 }
210 }
211 pumi_mesh->end(itr);
212
213 // Assign attributes for elements.
214 double ppt[3];
215 Vector cent(ppt, dim);
216 for (int el = 0; el < mesh->GetNE(); el++)
217 {
218 (mesh->GetElementTransformation(el))->
219 Transform(Geometries.GetCenter(mesh->GetElementBaseGeometry(el)),cent);
220 if (cent(0) <= -0.05)
221 {
222 mesh->SetAttribute(el, 1);
223 }
224 else if (cent(0) >= 0.05)
225 {
226 mesh->SetAttribute(el, 2);
227 }
228 else
229 {
230 mesh->SetAttribute(el, 3);
231 }
232 }
233 mesh->SetAttributes();
234 if (mesh->attributes.Max() < 2 || mesh->bdr_attributes.Max() < 2)
235 {
236 cerr << "\nInput mesh should have at least two materials and "
237 << "two boundary attributes! (See schematic in ex2.cpp)\n"
238 << endl;
239 return 3;
240 }
241
242 // 6. Refine the mesh to increase the resolution. In this example we do
243 // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
244 // largest number that gives a final mesh with no more than 5,000
245 // elements.
246 {
247 int ref_levels =
248 (int)floor(log(5000./mesh->GetNE())/log(2.)/dim);
249 for (int l = 0; l < ref_levels; l++)
250 {
251 mesh->UniformRefinement();
252 }
253 }
254
255 // 7. Define a finite element space on the mesh. Here we use vector finite
256 // elements, i.e. dim copies of a scalar finite element space. The vector
257 // dimension is specified by the last argument of the FiniteElementSpace
258 // constructor. For NURBS meshes, we use the (degree elevated) NURBS space
259 // associated with the mesh nodes.
261 FiniteElementSpace *fespace;
262 if (mesh->NURBSext)
263 {
264 fec = NULL;
265 fespace = mesh->GetNodes()->FESpace();
266 }
267 else
268 {
269 fec = new H1_FECollection(order, dim);
270 fespace = new FiniteElementSpace(mesh, fec, dim);
271 }
272 cout << "Number of finite element unknowns: " << fespace->GetTrueVSize()
273 << endl << "Assembling: " << flush;
274
275 // 8. Determine the list of true (i.e. conforming) essential boundary dofs.
276 // In this example, the boundary conditions are defined by marking only
277 // boundary attribute 1 from the mesh as essential and converting it to a
278 // list of true dofs.
279 Array<int> ess_tdof_list, ess_bdr(mesh->bdr_attributes.Max());
280 ess_bdr = 0;
281 ess_bdr[0] = 1;
282 fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
283
284 // 9. Set up the linear form b(.) which corresponds to the right-hand side of
285 // the FEM linear system. In this case, b_i equals the boundary integral
286 // of f*phi_i where f represents a "pull down" force on the Neumann part
287 // of the boundary and phi_i are the basis functions in the finite element
288 // fespace. The force is defined by the VectorArrayCoefficient object f,
289 // which is a vector of Coefficient objects. The fact that f is non-zero
290 // on boundary attribute 2 is indicated by the use of piece-wise constants
291 // coefficient for its last component.
293 for (int i = 0; i < dim-1; i++)
294 {
295 f.Set(i, new ConstantCoefficient(0.0));
296 }
297 {
298 Vector pull_force(mesh->bdr_attributes.Max());
299 pull_force = 0.0;
300 pull_force(1) = -3.0e-2;
301 f.Set(dim-1, new PWConstCoefficient(pull_force));
302 f.Set(dim-2, new PWConstCoefficient(pull_force));
303 }
304
305 LinearForm *b = new LinearForm(fespace);
306 b->AddBoundaryIntegrator(new VectorBoundaryLFIntegrator(f));
307 cout << "r.h.s. ... " << flush;
308 b->Assemble();
309
310 // 10. Define the solution vector x as a finite element grid function
311 // corresponding to fespace. Initialize x with initial guess of zero,
312 // which satisfies the boundary conditions.
313 GridFunction x(fespace);
314 x = 0.0;
315
316 // 11. Set up the bilinear form a(.,.) on the finite element space
317 // corresponding to the linear elasticity integrator with piece-wise
318 // constants coefficient lambda and mu.
319 Vector lambda(mesh->attributes.Max());
320 lambda = 1.0;
321 lambda(0) = lambda(1)*10;
322 lambda(1) = lambda(1)*100;
323 PWConstCoefficient lambda_func(lambda);
324 Vector mu(mesh->attributes.Max());
325 mu = 1.0;
326 mu(0) = mu(1)*10;
327 mu(1) = mu(1)*100;
328 PWConstCoefficient mu_func(mu);
329
330 BilinearForm *a = new BilinearForm(fespace);
331 a->AddDomainIntegrator(new ElasticityIntegrator(lambda_func,mu_func));
332
333 // 12. Assemble the bilinear form and the corresponding linear system,
334 // applying any necessary transformations such as: eliminating boundary
335 // conditions, applying conforming constraints for non-conforming AMR,
336 // static condensation, etc.
337 cout << "matrix ... " << flush;
338 if (static_cond) { a->EnableStaticCondensation(); }
339 a->Assemble();
340
341 SparseMatrix A;
342 Vector B, X;
343 a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
344 cout << "done." << endl;
345
346 cout << "Size of linear system: " << A.Height() << endl;
347
348#ifndef MFEM_USE_SUITESPARSE
349 // 13. Define a simple symmetric Gauss-Seidel preconditioner and use it to
350 // solve the system Ax=b with PCG.
351 GSSmoother M(A);
352 PCG(A, M, B, X, 1, 500, 1e-8, 0.0);
353#else
354 // 13. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
355 UMFPackSolver umf_solver;
356 umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
357 umf_solver.SetOperator(A);
358 umf_solver.Mult(B, X);
359#endif
360
361 // 14. Recover the solution as a finite element grid function.
362 a->RecoverFEMSolution(X, *b, x);
363
364 // 15. For non-NURBS meshes, make the mesh curved based on the finite element
365 // space. This means that we define the mesh elements through a fespace
366 // based transformation of the reference element. This allows us to save
367 // the displaced mesh as a curved mesh when using high-order finite
368 // element displacement field. We assume that the initial mesh (read from
369 // the file) is not higher order curved mesh compared to the chosen FE
370 // space.
371 if (!mesh->NURBSext)
372 {
373 mesh->SetNodalFESpace(fespace);
374 }
375
376 // 16. Save the displaced mesh and the inverted solution (which gives the
377 // backward displacements to the original grid). This output can be
378 // viewed later using GLVis: "glvis -m displaced.mesh -g sol.gf".
379 {
380 GridFunction *nodes = mesh->GetNodes();
381 *nodes += x;
382 x *= -1;
383 ofstream mesh_ofs("displaced.mesh");
384 mesh_ofs.precision(8);
385 mesh->Print(mesh_ofs);
386 ofstream sol_ofs("sol.gf");
387 sol_ofs.precision(8);
388 x.Save(sol_ofs);
389 }
390
391 // 17. Send the above data by socket to a GLVis server. Use the "n" and "b"
392 // keys in GLVis to visualize the displacements.
393 if (visualization)
394 {
395 char vishost[] = "localhost";
396 int visport = 19916;
397 socketstream sol_sock(vishost, visport);
398 sol_sock.precision(8);
399 sol_sock << "solution\n" << *mesh << x << flush;
400 }
401
402 // 18. Free the used memory.
403 delete a;
404 delete b;
405 if (fec)
406 {
407 delete fespace;
408 delete fec;
409 }
410 delete mesh;
411
412 pumi_mesh->destroyNative();
413 apf::destroyMesh(pumi_mesh);
414 PCU_Comm_Free();
415#ifdef MFEM_USE_SIMMETRIX
416 gmi_sim_stop();
417 Sim_unregisterAllKeys();
418#endif
419
420 return 0;
421}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:758
int Find(const T &el) const
Return the first index where 'el' is found; return -1 if not found.
Definition array.hpp:889
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
Definition array.cpp:23
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.
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
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.
const IntegrationPoint & GetCenter(int GeomType) const
Return the center of the given Geometry::Type, GeomType.
Definition geom.hpp:75
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
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:298
void SetAttribute(int i, int attr)
Set the attribute of element i.
Definition mesh.cpp:7629
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
const Element * GetBdrElement(int i) const
Return pointer to the i'th boundary element object.
Definition mesh.hpp:1354
void GetElementTransformation(int i, IsoparametricTransformation *ElTr) const
Builds the transformation defining the i-th element in ElTr. ElTr must be allocated in advance and wi...
Definition mesh.cpp:357
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:9294
virtual void SetNodalFESpace(FiniteElementSpace *nfes)
Definition mesh.cpp:6426
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:11295
Geometry::Type GetElementBaseGeometry(int i) const
Definition mesh.hpp:1455
Array< int > attributes
A list of all unique element attributes used by the Mesh.
Definition mesh.hpp:288
virtual void SetAttributes()
Determine the sets of unique attribute values in domain and boundary elements.
Definition mesh.cpp:1819
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.
A piecewise constant coefficient with the constants keyed off the element attribute numbers.
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 coefficient defined by an array of scalar coefficients. Coefficients that are not set will eva...
Vector data type.
Definition vector.hpp:82
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:337
int dim
Definition ex24.cpp:53
real_t mu
Definition ex25.cpp:140
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
Geometry Geometries
Definition fe.cpp:49
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
void filter_dos(std::string &line)
Check for, and remove, a trailing '\r' from and std::string.
Definition text.hpp:45
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void skip_comment_lines(std::istream &is, const char comment_char)
Check if the stream starts with comment_char. If so skip it.
Definition text.hpp:31
const char vishost[]
std::array< int, NCMesh::MaxFaceNodes > nodes