MFEM  v4.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 - 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 // ex1 -m ../../data/beam-hex-nurbs.mesh -pm 1 -ps 2
23 //
24 // Description: This example code demonstrates the use of MFEM to define a
25 // simple finite element discretization of the Laplace problem
26 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
27 // Specifically, we discretize using a FE space of the specified
28 // order, or if order < 1 using an isoparametric/isogeometric
29 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
30 // NURBS mesh, etc.)
31 //
32 // The example highlights the use of mesh refinement, finite
33 // element grid functions, as well as linear and bilinear forms
34 // corresponding to the left-hand side and right-hand side of the
35 // discrete linear system. We also cover the explicit elimination
36 // of essential boundary conditions, static condensation, and the
37 // optional connection to the GLVis tool for visualization.
38 
39 #include "mfem.hpp"
40 #include <fstream>
41 #include <iostream>
42 
43 using namespace std;
44 using namespace mfem;
45 
46 int main(int argc, char *argv[])
47 {
48  // 1. Parse command-line options.
49  const char *mesh_file = "../../data/star.mesh";
50  const char *per_file = "none";
51  Array<int> master(0);
52  Array<int> slave(0);
53  bool static_cond = false;
54  bool visualization = 1;
55  Array<int> order(1);
56  order[0] = 1;
57 
58  OptionsParser args(argc, argv);
59  args.AddOption(&mesh_file, "-m", "--mesh",
60  "Mesh file to use.");
61  args.AddOption(&per_file, "-p", "--per",
62  "Periodic BCS file.");
63  args.AddOption(&master, "-pm", "--master",
64  "Master boundaries for periodic BCs");
65  args.AddOption(&slave, "-ps", "--slave",
66  "Slave boundaries for periodic BCs");
67  args.AddOption(&order, "-o", "--order",
68  "Finite element order (polynomial degree) or -1 for"
69  " isoparametric space.");
70  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
71  "--no-static-condensation", "Enable static condensation.");
72  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
73  "--no-visualization",
74  "Enable or disable GLVis visualization.");
75  args.Parse();
76  if (!args.Good())
77  {
78  args.PrintUsage(cout);
79  return 1;
80  }
81  args.PrintOptions(cout);
82 
83  // 2. Read the mesh from the given mesh file. We can handle triangular,
84  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
85  // the same code.
86  Mesh *mesh = new Mesh(mesh_file, 1, 1);
87  int dim = mesh->Dimension();
88 
89  // 3. Refine the mesh to increase the resolution. In this example we do
90  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
91  // largest number that gives a final mesh with no more than 50,000
92  // elements.
93  {
94  int ref_levels =
95  (int)floor(log(50000./mesh->GetNE())/log(2.)/dim);
96  for (int l = 0; l < ref_levels; l++)
97  {
98  mesh->UniformRefinement();
99  }
100  }
101 
102  // 4. Define a finite element space on the mesh. Here we use continuous
103  // Lagrange finite elements of the specified order. If order < 1, we
104  // instead use an isoparametric/isogeometric space.
106  NURBSExtension *NURBSext = NULL;
107  int own_fec = 0;
108 
109  if (mesh->NURBSext)
110  {
111  fec = new NURBSFECollection(order[0]);
112  own_fec = 1;
113 
114  int nkv = mesh->NURBSext->GetNKV();
115  if (order.Size() == 1)
116  {
117  int tmp = order[0];
118  order.SetSize(nkv);
119  order = tmp;
120  }
121 
122  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
123  NURBSext = new NURBSExtension(mesh->NURBSext, order);
124 
125  // Read periodic BCs from file
126  std::ifstream in;
127  in.open(per_file, std::ifstream::in);
128  if (in.is_open())
129  {
130  int psize;
131  in >> psize;
132  master.SetSize(psize);
133  slave.SetSize(psize);
134  master.Load(in, psize);
135  slave.Load(in, psize);
136  in.close();
137  }
138  master.Print();
139  slave.Print();
140  NURBSext->ConnectBoundaries(master,slave);
141  }
142  else if (order[0] == -1) // Isoparametric
143  {
144  if (mesh->GetNodes())
145  {
146  fec = mesh->GetNodes()->OwnFEC();
147  own_fec = 0;
148  cout << "Using isoparametric FEs: " << fec->Name() << endl;
149  }
150  else
151  {
152  cout <<"Mesh does not have FEs --> Assume order 1.\n";
153  fec = new H1_FECollection(1, dim);
154  own_fec = 1;
155  }
156  }
157  else
158  {
159  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
160  fec = new H1_FECollection(abs(order[0]), dim);
161  own_fec = 1;
162  }
163  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, NURBSext, fec);
164  cout << "Number of finite element unknowns: "
165  << fespace->GetTrueVSize() << endl;
166 
167  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
168  // In this example, the boundary conditions are defined by marking all
169  // the boundary attributes from the mesh as essential (Dirichlet) and
170  // converting them to a list of true dofs.
171  Array<int> ess_tdof_list;
172  if (mesh->bdr_attributes.Size())
173  {
174  Array<int> ess_bdr(mesh->bdr_attributes.Max());
175  ess_bdr = 1;
176  // Remove periodic BCs
177  for (int i = 0; i < master.Size(); i++)
178  {
179  ess_bdr[master[i]-1] = 0;
180  ess_bdr[slave[i]-1] = 0;
181  }
182  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
183  }
184 
185  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
186  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
187  // the basis functions in the finite element fespace.
188  LinearForm *b = new LinearForm(fespace);
189  ConstantCoefficient one(1.0);
191  b->Assemble();
192 
193  // 7. Define the solution vector x as a finite element grid function
194  // corresponding to fespace. Initialize x with initial guess of zero,
195  // which satisfies the boundary conditions.
196  GridFunction x(fespace);
197  x = 0.0;
198 
199  // 8. Set up the bilinear form a(.,.) on the finite element space
200  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
201  // domain integrator.
202  BilinearForm *a = new BilinearForm(fespace);
204 
205  // 9. Assemble the bilinear form and the corresponding linear system,
206  // applying any necessary transformations such as: eliminating boundary
207  // conditions, applying conforming constraints for non-conforming AMR,
208  // static condensation, etc.
209  if (static_cond) { a->EnableStaticCondensation(); }
210  a->Assemble();
211 
212  SparseMatrix A;
213  Vector B, X;
214  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
215 
216  cout << "Size of linear system: " << A.Height() << endl;
217 
218 #ifndef MFEM_USE_SUITESPARSE
219  // 10. Define a simple symmetric Gauss-Seidel preconditioner and use it to
220  // solve the system A X = B with PCG.
221  GSSmoother M(A);
222  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
223 #else
224  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
225  UMFPackSolver umf_solver;
226  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
227  umf_solver.SetOperator(A);
228  umf_solver.Mult(B, X);
229 #endif
230 
231  // 11. Recover the solution as a finite element grid function.
232  a->RecoverFEMSolution(X, *b, x);
233 
234  // 12. Save the refined mesh and the solution. This output can be viewed later
235  // using GLVis: "glvis -m refined.mesh -g sol.gf".
236  ofstream mesh_ofs("refined.mesh");
237  mesh_ofs.precision(8);
238  mesh->Print(mesh_ofs);
239  ofstream sol_ofs("sol.gf");
240  sol_ofs.precision(8);
241  x.Save(sol_ofs);
242 
243  // 13. Send the solution by socket to a GLVis server.
244  if (visualization)
245  {
246  char vishost[] = "localhost";
247  int visport = 19916;
248  socketstream sol_sock(vishost, visport);
249  sol_sock.precision(8);
250  sol_sock << "solution\n" << *mesh << x << flush;
251  }
252 
253  // 14. Save data in the VisIt format
254  VisItDataCollection visit_dc("Example1", mesh);
255  visit_dc.RegisterField("solution", &x);
256  visit_dc.Save();
257 
258  // 15. Free the used memory.
259  delete a;
260  delete b;
261  delete fespace;
262  if (own_fec) { delete fec; }
263  delete mesh;
264 
265  return 0;
266 }
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:278
int Size() const
Logical size of the array.
Definition: array.hpp:118
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1156
void Load(std::istream &in, int fmt=0)
Read an Array from the stream in using format fmt. The format fmt can be:
Definition: array.cpp:53
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(...
int GetNKV() const
Definition: nurbs.hpp:342
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:676
virtual void Save()
Save the collection and a VisIt root file.
int main(int argc, char *argv[])
Definition: ex1.cpp:58
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:2627
int dim
Definition: ex3.cpp:48
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: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:7591
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: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:461
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:350
int Dimension() const
Definition: mesh.hpp:713
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:179
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:85
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
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:618
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:181
virtual const char * Name() const
Definition: fe_coll.hpp:53
void ConnectBoundaries()
Definition: nurbs.cpp:1640
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void Print(std::ostream &out=mfem::out, int width=4) const
Prints array to stream with width elements per row.
Definition: array.cpp:23
Vector data type.
Definition: vector.hpp:48
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:5837
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:1794
void EnableStaticCondensation()
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:1697
bool Good() const
Definition: optparser.hpp:122