MFEM  v3.4
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 - NURBS Version
2 //
3 // Compile with: make ex1
4 //
5 // Sample runs: ex1 -m ../../data/square-disc.mesh
6 // ex1 -m ../../data/star.mesh
7 // ex1 -m ../../data/escher.mesh
8 // ex1 -m ../../data/fichera.mesh
9 // ex1 -m ../../data/square-disc-p2.vtk -o 2
10 // ex1 -m ../../data/square-disc-p3.mesh -o 3
11 // ex1 -m ../../data/square-disc-nurbs.mesh -o -1
12 // ex1 -m ../../data/disc-nurbs.mesh -o -1
13 // ex1 -m ../../data/pipe-nurbs.mesh -o -1
14 // ex1 -m ../../data/star-surf.mesh
15 // ex1 -m ../../data/square-disc-surf.mesh
16 // ex1 -m ../../data/inline-segment.mesh
17 // ex1 -m ../../data/amr-quad.mesh
18 // ex1 -m ../../data/amr-hex.mesh
19 // ex1 -m ../../data/fichera-amr.mesh
20 // ex1 -m ../../data/mobius-strip.mesh
21 // ex1 -m ../../data/mobius-strip.mesh -o -1 -sc
22 //
23 // Description: This example code demonstrates the use of MFEM to define a
24 // simple finite element discretization of the Laplace problem
25 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
26 // Specifically, we discretize using a FE space of the specified
27 // order, or if order < 1 using an isoparametric/isogeometric
28 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
29 // NURBS mesh, etc.)
30 //
31 // The example highlights the use of mesh refinement, finite
32 // element grid functions, as well as linear and bilinear forms
33 // corresponding to the left-hand side and right-hand side of the
34 // discrete linear system. We also cover the explicit elimination
35 // of essential boundary conditions, static condensation, and the
36 // optional connection to the GLVis tool for visualization.
37 
38 #include "mfem.hpp"
39 #include <fstream>
40 #include <iostream>
41 
42 using namespace std;
43 using namespace mfem;
44 
45 int main(int argc, char *argv[])
46 {
47  // 1. Parse command-line options.
48  const char *mesh_file = "../../data/star.mesh";
49  bool static_cond = false;
50  bool visualization = 1;
51  Array<int> order(1);
52  order[0] = 1;
53 
54  OptionsParser args(argc, argv);
55  args.AddOption(&mesh_file, "-m", "--mesh",
56  "Mesh file to use.");
57  args.AddOption(&order, "-o", "--order",
58  "Finite element order (polynomial degree) or -1 for"
59  " isoparametric space.");
60  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
61  "--no-static-condensation", "Enable static condensation.");
62  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
63  "--no-visualization",
64  "Enable or disable GLVis visualization.");
65  args.Parse();
66  if (!args.Good())
67  {
68  args.PrintUsage(cout);
69  return 1;
70  }
71  args.PrintOptions(cout);
72 
73  // 2. Read the mesh from the given mesh file. We can handle triangular,
74  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
75  // the same code.
76  Mesh *mesh = new Mesh(mesh_file, 1, 1);
77  int dim = mesh->Dimension();
78 
79  // 3. Refine the mesh to increase the resolution. In this example we do
80  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
81  // largest number that gives a final mesh with no more than 50,000
82  // elements.
83  {
84  int ref_levels =
85  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
86  for (int l = 0; l < ref_levels; l++)
87  {
88  mesh->UniformRefinement();
89  }
90  }
91 
92  // 4. Define a finite element space on the mesh. Here we use continuous
93  // Lagrange finite elements of the specified order. If order < 1, we
94  // instead use an isoparametric/isogeometric space.
96  NURBSExtension *NURBSext = NULL;
97  int own_fec = 0;
98 
99  if (order[0] == -1) // Isoparametric
100  {
101  if (mesh->GetNodes())
102  {
103  fec = mesh->GetNodes()->OwnFEC();
104  own_fec = 0;
105  cout << "Using isoparametric FEs: " << fec->Name() << endl;
106  }
107  else
108  {
109  cout <<"Mesh does not have FEs --> Assume order 1.\n";
110  fec = new H1_FECollection(1, dim);
111  own_fec = 1;
112  }
113  }
114  else if (mesh->NURBSext && (order[0] > 0) ) // Subparametric NURBS
115  {
116  fec = new NURBSFECollection(order[0]);
117  own_fec = 1;
118  int nkv = mesh->NURBSext->GetNKV();
119 
120  if (order.Size() == 1)
121  {
122  int tmp = order[0];
123  order.SetSize(nkv);
124  order = tmp;
125  }
126  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
127  NURBSext = new NURBSExtension(mesh->NURBSext, order);
128  }
129  else
130  {
131  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
132  fec = new H1_FECollection(abs(order[0]), dim);
133  own_fec = 1;
134  }
135  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, NURBSext, fec);
136  cout << "Number of finite element unknowns: "
137  << fespace->GetTrueVSize() << endl;
138 
139  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
140  // In this example, the boundary conditions are defined by marking all
141  // the boundary attributes from the mesh as essential (Dirichlet) and
142  // converting them to a list of true dofs.
143  Array<int> ess_tdof_list;
144  if (mesh->bdr_attributes.Size())
145  {
146  Array<int> ess_bdr(mesh->bdr_attributes.Max());
147  ess_bdr = 1;
148  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
149  }
150 
151  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
152  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
153  // the basis functions in the finite element fespace.
154  LinearForm *b = new LinearForm(fespace);
155  ConstantCoefficient one(1.0);
157  b->Assemble();
158 
159  // 7. Define the solution vector x as a finite element grid function
160  // corresponding to fespace. Initialize x with initial guess of zero,
161  // which satisfies the boundary conditions.
162  GridFunction x(fespace);
163  x = 0.0;
164 
165  // 8. Set up the bilinear form a(.,.) on the finite element space
166  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
167  // domain integrator.
168  BilinearForm *a = new BilinearForm(fespace);
170 
171  // 9. Assemble the bilinear form and the corresponding linear system,
172  // applying any necessary transformations such as: eliminating boundary
173  // conditions, applying conforming constraints for non-conforming AMR,
174  // static condensation, etc.
175  if (static_cond) { a->EnableStaticCondensation(); }
176  a->Assemble();
177 
178  SparseMatrix A;
179  Vector B, X;
180  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
181 
182  cout << "Size of linear system: " << A.Height() << endl;
183 
184 #ifndef MFEM_USE_SUITESPARSE
185  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
186  // solve the system A X = B with PCG.
187  GSSmoother M(A);
188  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
189 #else
190  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
191  UMFPackSolver umf_solver;
192  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
193  umf_solver.SetOperator(A);
194  umf_solver.Mult(B, X);
195 #endif
196 
197  // 11. Recover the solution as a finite element grid function.
198  a->RecoverFEMSolution(X, *b, x);
199 
200  // 12. Save the refined mesh and the solution. This output can be viewed later
201  // using GLVis: "glvis -m refined.mesh -g sol.gf".
202  ofstream mesh_ofs("refined.mesh");
203  mesh_ofs.precision(8);
204  mesh->Print(mesh_ofs);
205  ofstream sol_ofs("sol.gf");
206  sol_ofs.precision(8);
207  x.Save(sol_ofs);
208 
209  // 13. Send the solution by socket to a GLVis server.
210  if (visualization)
211  {
212  char vishost[] = "localhost";
213  int visport = 19916;
214  socketstream sol_sock(vishost, visport);
215  sol_sock.precision(8);
216  sol_sock << "solution\n" << *mesh << x << flush;
217  }
218 
219  // 14. Save data in the VisIt format
220  VisItDataCollection visit_dc("Example1", mesh);
221  visit_dc.RegisterField("solution", &x);
222  visit_dc.Save();
223 
224  // 15. Free the used memory.
225  delete a;
226  delete b;
227  delete fespace;
228  if (own_fec) { delete fec; }
229  delete mesh;
230 
231  return 0;
232 }
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:266
int Size() const
Logical size of the array.
Definition: array.hpp:133
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1034
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:57
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:51
int GetNKV() const
Definition: nurbs.hpp:319
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:365
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
virtual void Save()
Save the collection and a VisIt root file.
int main(int argc, char *argv[])
Definition: ex1.cpp:45
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:344
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2349
int dim
Definition: ex3.cpp:47
Data type sparse matrix.
Definition: sparsemat.hpp:38
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:36
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:146
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
Data collection with VisIt I/O routines.
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:109
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:457
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:256
int Dimension() const
Definition: mesh.hpp:645
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator.
Definition: linearform.cpp:19
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:174
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:355
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:66
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:74
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:569
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:176
virtual const char * Name() const
Definition: fe_coll.hpp:50
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator.
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:5422
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:79
Class for linear form - Vector with associated FE space and LFIntegrators.
Definition: linearform.hpp:23
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, SparseMatrix &A, Vector &X, Vector &B, int copy_interior=0)
Form a linear system, A X = B.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1785
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1688
bool Good() const
Definition: optparser.hpp:120