MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nurbs_ex1.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - NURBS Version
2 //
3 // Compile with: make nurbs_ex1
4 //
5 // Sample runs: nurbs_ex1 -m square-nurbs.mesh -o 2 -no-ibp
6 // nurbs_ex1 -m cube-nurbs.mesh -o 2 -no-ibp
7 // nurbs_ex1 -m pipe-nurbs-2d.mesh -o 2 -no-ibp
8 // nurbs_ex1 -m ../../data/square-disc-nurbs.mesh -o -1
9 // nurbs_ex1 -m ../../data/disc-nurbs.mesh -o -1
10 // nurbs_ex1 -m ../../data/pipe-nurbs.mesh -o -1
11 // nurbs_ex1 -m ../../data/beam-hex-nurbs.mesh -pm 1 -ps 2
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 #include "mfem.hpp"
29 #include <fstream>
30 #include <iostream>
31 
32 using namespace std;
33 using namespace mfem;
34 
35 
36 /** Class for integrating the bilinear form a(u,v) := (Q Laplace u, v) where Q
37  can be a scalar coefficient. */
38 class Diffusion2Integrator: public BilinearFormIntegrator
39 {
40 private:
41 #ifndef MFEM_THREAD_SAFE
42  Vector shape,laplace;
43 #endif
44  Coefficient *Q;
45 
46 public:
47  /// Construct a diffusion integrator with coefficient Q = 1
48  Diffusion2Integrator() { Q = NULL; }
49 
50  /// Construct a diffusion integrator with a scalar coefficient q
51  Diffusion2Integrator (Coefficient &q) : Q(&q) { }
52 
53  /** Given a particular Finite Element
54  computes the element stiffness matrix elmat. */
55  virtual void AssembleElementMatrix(const FiniteElement &el,
57  DenseMatrix &elmat)
58  {
59  int nd = el.GetDof();
60  int dim = el.GetDim();
61  double w;
62 
63 #ifdef MFEM_THREAD_SAFE
64  Vector shape[nd];
65  Vector laplace(nd);
66 #else
67  shape.SetSize(nd);
68  laplace.SetSize(nd);
69 #endif
70  elmat.SetSize(nd);
71 
72  const IntegrationRule *ir = IntRule;
73  if (ir == NULL)
74  {
75  int order;
76  if (el.Space() == FunctionSpace::Pk)
77  {
78  order = 2*el.GetOrder() - 2;
79  }
80  else
81  {
82  order = 2*el.GetOrder() + dim - 1;
83  }
84 
85  if (el.Space() == FunctionSpace::rQk)
86  {
87  ir = &RefinedIntRules.Get(el.GetGeomType(),order);
88  }
89  else
90  {
91  ir = &IntRules.Get(el.GetGeomType(),order);
92  }
93  }
94 
95  elmat = 0.0;
96  for (int i = 0; i < ir->GetNPoints(); i++)
97  {
98  const IntegrationPoint &ip = ir->IntPoint(i);
99  Trans.SetIntPoint(&ip);
100  w = -ip.weight * Trans.Weight();
101 
102  el.CalcShape(ip, shape);
103  el.CalcPhysLaplacian(Trans, laplace);
104 
105  if (Q)
106  {
107  w *= Q->Eval(Trans, ip);
108  }
109 
110  for (int j = 0; j < nd; j++)
111  {
112  for (int i = 0; i < nd; i++)
113  {
114  elmat(i, j) += w*shape(i)*laplace(j);
115  }
116  }
117  }
118  }
119 
120 };
121 
122 int main(int argc, char *argv[])
123 {
124  // 1. Parse command-line options.
125  const char *mesh_file = "../../data/star.mesh";
126  const char *per_file = "none";
127  Array<int> master(0);
128  Array<int> slave(0);
129  bool static_cond = false;
130  bool visualization = 1;
131  bool ibp = 1;
132  Array<int> order(1);
133  order[0] = 1;
134 
135  OptionsParser args(argc, argv);
136  args.AddOption(&mesh_file, "-m", "--mesh",
137  "Mesh file to use.");
138  args.AddOption(&per_file, "-p", "--per",
139  "Periodic BCS file.");
140  args.AddOption(&master, "-pm", "--master",
141  "Master boundaries for periodic BCs");
142  args.AddOption(&slave, "-ps", "--slave",
143  "Slave boundaries for periodic BCs");
144  args.AddOption(&order, "-o", "--order",
145  "Finite element order (polynomial degree) or -1 for"
146  " isoparametric space.");
147  args.AddOption(&ibp, "-ibp", "--ibp", "-no-ibp",
148  "--no-ibp",
149  "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
150  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
151  "--no-static-condensation", "Enable static condensation.");
152  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
153  "--no-visualization",
154  "Enable or disable GLVis visualization.");
155  args.Parse();
156  if (!args.Good())
157  {
158  args.PrintUsage(cout);
159  return 1;
160  }
161  args.PrintOptions(cout);
162 
163  // 2. Read the mesh from the given mesh file. We can handle triangular,
164  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
165  // the same code.
166  Mesh *mesh = new Mesh(mesh_file, 1, 1);
167  int dim = mesh->Dimension();
168 
169  // 3. Refine the mesh to increase the resolution. In this example we do
170  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
171  // largest number that gives a final mesh with no more than 50,000
172  // elements.
173  {
174  int ref_levels =
175  (int)floor(log(5000./mesh->GetNE())/log(2.)/dim);
176  for (int l = 0; l < ref_levels; l++)
177  {
178  mesh->UniformRefinement();
179  }
180  }
181 
182  // 4. Define a finite element space on the mesh. Here we use continuous
183  // Lagrange finite elements of the specified order. If order < 1, we
184  // instead use an isoparametric/isogeometric space.
186  NURBSExtension *NURBSext = NULL;
187  int own_fec = 0;
188 
189  if (mesh->NURBSext)
190  {
191  fec = new NURBSFECollection(order[0]);
192  own_fec = 1;
193 
194  int nkv = mesh->NURBSext->GetNKV();
195  if (order.Size() == 1)
196  {
197  int tmp = order[0];
198  order.SetSize(nkv);
199  order = tmp;
200  }
201 
202  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
203  NURBSext = new NURBSExtension(mesh->NURBSext, order);
204 
205  // Read periodic BCs from file
206  std::ifstream in;
207  in.open(per_file, std::ifstream::in);
208  if (in.is_open())
209  {
210  int psize;
211  in >> psize;
212  master.SetSize(psize);
213  slave.SetSize(psize);
214  master.Load(in, psize);
215  slave.Load(in, psize);
216  in.close();
217  }
218  master.Print();
219  slave.Print();
220  NURBSext->ConnectBoundaries(master,slave);
221  }
222  else if (order[0] == -1) // Isoparametric
223  {
224  if (mesh->GetNodes())
225  {
226  fec = mesh->GetNodes()->OwnFEC();
227  own_fec = 0;
228  cout << "Using isoparametric FEs: " << fec->Name() << endl;
229  }
230  else
231  {
232  cout <<"Mesh does not have FEs --> Assume order 1.\n";
233  fec = new H1_FECollection(1, dim);
234  own_fec = 1;
235  }
236  }
237  else
238  {
239  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
240  fec = new H1_FECollection(abs(order[0]), dim);
241  own_fec = 1;
242  }
243 
244  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, NURBSext, fec);
245  cout << "Number of finite element unknowns: "
246  << fespace->GetTrueVSize() << endl;
247 
248  if (!ibp)
249  {
250  if (!mesh->NURBSext)
251  {
252  cout << "No integration by parts requires a NURBS mesh."<< endl;
253  return 2;
254  }
255  if (mesh->NURBSext->GetNP()>1)
256  {
257  cout << "No integration by parts requires a NURBS mesh, with only 1 patch."<<
258  endl;
259  cout << "A C_1 discretisation is required."<< endl;
260  cout << "Currently only C_0 multipatch coupling implemented."<< endl;
261  return 3;
262  }
263  if (order[0]<2)
264  {
265  cout << "No integration by parts requires at least quadratic NURBS."<< endl;
266  cout << "A C_1 discretisation is required."<< endl;
267  return 4;
268  }
269  }
270 
271 
272 
273  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
274  // In this example, the boundary conditions are defined by marking all
275  // the boundary attributes from the mesh as essential (Dirichlet) and
276  // converting them to a list of true dofs.
277  Array<int> ess_tdof_list;
278  if (mesh->bdr_attributes.Size())
279  {
280  Array<int> ess_bdr(mesh->bdr_attributes.Max());
281  ess_bdr = 1;
282  // Remove periodic BCs
283  for (int i = 0; i < master.Size(); i++)
284  {
285  ess_bdr[master[i]-1] = 0;
286  ess_bdr[slave[i]-1] = 0;
287  }
288  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
289  }
290 
291  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
292  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
293  // the basis functions in the finite element fespace.
294  LinearForm *b = new LinearForm(fespace);
295  ConstantCoefficient one(1.0);
297  b->Assemble();
298 
299  // 7. Define the solution vector x as a finite element grid function
300  // corresponding to fespace. Initialize x with initial guess of zero,
301  // which satisfies the boundary conditions.
302  GridFunction x(fespace);
303  x = 0.0;
304 
305  // 8. Set up the bilinear form a(.,.) on the finite element space
306  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
307  // domain integrator.
308  BilinearForm *a = new BilinearForm(fespace);
309  if (ibp)
310  {
312  }
313  else
314  {
315  a->AddDomainIntegrator(new Diffusion2Integrator(one));
316  }
317 
318  // 9. Assemble the bilinear form and the corresponding linear system,
319  // applying any necessary transformations such as: eliminating boundary
320  // conditions, applying conforming constraints for non-conforming AMR,
321  // static condensation, etc.
322  if (static_cond) { a->EnableStaticCondensation(); }
323  a->Assemble();
324 
325  SparseMatrix A;
326  Vector B, X;
327  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
328 
329  cout << "Size of linear system: " << A.Height() << endl;
330 
331 #ifndef MFEM_USE_SUITESPARSE
332  // 10. Define a simple Jacobi preconditioner and use it to
333  // solve the system A X = B with PCG.
334  GSSmoother M(A);
335  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
336 #else
337  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
338  UMFPackSolver umf_solver;
339  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
340  umf_solver.SetOperator(A);
341  umf_solver.Mult(B, X);
342 #endif
343 
344  // 11. Recover the solution as a finite element grid function.
345  a->RecoverFEMSolution(X, *b, x);
346 
347  // 12. Save the refined mesh and the solution. This output can be viewed later
348  // using GLVis: "glvis -m refined.mesh -g sol.gf".
349  ofstream mesh_ofs("refined.mesh");
350  mesh_ofs.precision(8);
351  mesh->Print(mesh_ofs);
352  ofstream sol_ofs("sol.gf");
353  sol_ofs.precision(8);
354  x.Save(sol_ofs);
355 
356  // 13. Send the solution by socket to a GLVis server.
357  if (visualization)
358  {
359  char vishost[] = "localhost";
360  int visport = 19916;
361  socketstream sol_sock(vishost, visport);
362  sol_sock.precision(8);
363  sol_sock << "solution\n" << *mesh << x << flush;
364  }
365 
366  // 14. Save data in the VisIt format
367  VisItDataCollection visit_dc("Example1", mesh);
368  visit_dc.RegisterField("solution", &x);
369  visit_dc.Save();
370 
371  // 15. Free the used memory.
372  delete a;
373  delete b;
374  delete fespace;
375  if (own_fec) { delete fec; }
376  delete mesh;
377 
378  return 0;
379 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
Abstract class for Finite Elements.
Definition: fe.hpp:232
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:289
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
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:308
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 an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:890
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 SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
virtual void CalcPhysLaplacian(ElementTransformation &Trans, Vector &Laplacian) const
Evaluate the Laplacian of all shape functions of a scalar finite element in reference space at the gi...
Definition: fe.cpp:206
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:56
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:350
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
Recover the solution of a linear system formed with FormLinearSystem().
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:318
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
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
virtual void Save()
Save the collection and a VisIt root file.
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:328
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
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:311
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
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
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:153
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
Data collection with VisIt I/O routines.
virtual void CalcShape(const IntegrationPoint &ip, Vector &shape) const =0
Evaluate the values of all shape functions of a scalar finite element in reference space at the given...
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
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:24
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
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:314
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
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
int GetNP() const
Definition: nurbs.hpp:341
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:191
virtual const char * Name() const
Definition: fe_coll.hpp:53
Class for integration point with weight.
Definition: intrules.hpp:25
void ConnectBoundaries()
Definition: nurbs.cpp:1793
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
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
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition: intrules.hpp:379
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
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:90
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:376
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