MFEM  v4.4.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 ../../data/square-nurbs.mesh -o 2 -no-ibp
6 // nurbs_ex1 -m ../../data/square-nurbs.mesh -o 2 --weak-bc
7 // nurbs_ex1 -m ../../data/cube-nurbs.mesh -o 2 -no-ibp
8 // nurbs_ex1 -m ../../data/pipe-nurbs-2d.mesh -o 2 -no-ibp
9 // nurbs_ex1 -m ../../data/square-disc-nurbs.mesh -o -1
10 // nurbs_ex1 -m ../../data/disc-nurbs.mesh -o -1
11 // nurbs_ex1 -m ../../data/pipe-nurbs.mesh -o -1
12 // nurbs_ex1 -m ../../data/beam-hex-nurbs.mesh -pm 1 -ps 2
13 //
14 // Description: This example code demonstrates the use of MFEM to define a
15 // simple finite element discretization of the Laplace problem
16 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
17 // The boundary conditions can be enforced either strongly or weakly.
18 // Specifically, we discretize using a FE space of the specified
19 // order, or if order < 1 using an isoparametric/isogeometric
20 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
21 // NURBS mesh, etc.)
22 //
23 // The example highlights the use of mesh refinement, finite
24 // element grid functions, as well as linear and bilinear forms
25 // corresponding to the left-hand side and right-hand side of the
26 // discrete linear system. We also cover the explicit elimination
27 // of essential boundary conditions, static condensation, and the
28 // optional connection to the GLVis tool for visualization.
29 
30 #include "mfem.hpp"
31 #include <fstream>
32 #include <iostream>
33 
34 using namespace std;
35 using namespace mfem;
36 
37 
38 /** Class for integrating the bilinear form a(u,v) := (Q Laplace u, v) where Q
39  can be a scalar coefficient. */
40 class Diffusion2Integrator: public BilinearFormIntegrator
41 {
42 private:
43 #ifndef MFEM_THREAD_SAFE
44  Vector shape,laplace;
45 #endif
46  Coefficient *Q;
47 
48 public:
49  /// Construct a diffusion integrator with coefficient Q = 1
50  Diffusion2Integrator() { Q = NULL; }
51 
52  /// Construct a diffusion integrator with a scalar coefficient q
53  Diffusion2Integrator (Coefficient &q) : Q(&q) { }
54 
55  /** Given a particular Finite Element
56  computes the element stiffness matrix elmat. */
57  virtual void AssembleElementMatrix(const FiniteElement &el,
59  DenseMatrix &elmat)
60  {
61  int nd = el.GetDof();
62  int dim = el.GetDim();
63  double w;
64 
65 #ifdef MFEM_THREAD_SAFE
66  Vector shape(nd);
67  Vector laplace(nd);
68 #else
69  shape.SetSize(nd);
70  laplace.SetSize(nd);
71 #endif
72  elmat.SetSize(nd);
73 
74  const IntegrationRule *ir = IntRule;
75  if (ir == NULL)
76  {
77  int order;
78  if (el.Space() == FunctionSpace::Pk)
79  {
80  order = 2*el.GetOrder() - 2;
81  }
82  else
83  {
84  order = 2*el.GetOrder() + dim - 1;
85  }
86 
87  if (el.Space() == FunctionSpace::rQk)
88  {
89  ir = &RefinedIntRules.Get(el.GetGeomType(),order);
90  }
91  else
92  {
93  ir = &IntRules.Get(el.GetGeomType(),order);
94  }
95  }
96 
97  elmat = 0.0;
98  for (int i = 0; i < ir->GetNPoints(); i++)
99  {
100  const IntegrationPoint &ip = ir->IntPoint(i);
101  Trans.SetIntPoint(&ip);
102  w = -ip.weight * Trans.Weight();
103 
104  el.CalcShape(ip, shape);
105  el.CalcPhysLaplacian(Trans, laplace);
106 
107  if (Q)
108  {
109  w *= Q->Eval(Trans, ip);
110  }
111 
112  for (int jj = 0; jj < nd; jj++)
113  {
114  for (int ii = 0; ii < nd; ii++)
115  {
116  elmat(ii, jj) += w*shape(ii)*laplace(jj);
117  }
118  }
119  }
120  }
121 
122 };
123 
124 int main(int argc, char *argv[])
125 {
126  // 1. Parse command-line options.
127  const char *mesh_file = "../../data/star.mesh";
128  const char *per_file = "none";
129  int ref_levels = -1;
130  Array<int> master(0);
131  Array<int> slave(0);
132  bool static_cond = false;
133  bool visualization = 1;
134  bool ibp = 1;
135  bool strongBC = 1;
136  double kappa = -1;
137  Array<int> order(1);
138  order[0] = 1;
139 
140  OptionsParser args(argc, argv);
141  args.AddOption(&mesh_file, "-m", "--mesh",
142  "Mesh file to use.");
143  args.AddOption(&ref_levels, "-r", "--refine",
144  "Number of times to refine the mesh uniformly, -1 for auto.");
145  args.AddOption(&per_file, "-p", "--per",
146  "Periodic BCS file.");
147  args.AddOption(&master, "-pm", "--master",
148  "Master boundaries for periodic BCs");
149  args.AddOption(&slave, "-ps", "--slave",
150  "Slave boundaries for periodic BCs");
151  args.AddOption(&order, "-o", "--order",
152  "Finite element order (polynomial degree) or -1 for"
153  " isoparametric space.");
154  args.AddOption(&ibp, "-ibp", "--ibp", "-no-ibp",
155  "--no-ibp",
156  "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
157  args.AddOption(&strongBC, "-sbc", "--strong-bc", "-wbc",
158  "--weak-bc",
159  "Selects strong or weak enforcement of Dirichlet BCs.");
160  args.AddOption(&kappa, "-k", "--kappa",
161  "Sets the SIPG penalty parameters, should be positive."
162  " Negative values are replaced with (order+1)^2.");
163  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
164  "--no-static-condensation", "Enable static condensation.");
165  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
166  "--no-visualization",
167  "Enable or disable GLVis visualization.");
168  args.Parse();
169  if (!args.Good())
170  {
171  args.PrintUsage(cout);
172  return 1;
173  }
174  if (!strongBC & (kappa < 0))
175  {
176  kappa = 4*(order.Max()+1)*(order.Max()+1);
177  }
178  args.PrintOptions(cout);
179 
180  // 2. Read the mesh from the given mesh file. We can handle triangular,
181  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes with
182  // the same code.
183  Mesh *mesh = new Mesh(mesh_file, 1, 1);
184  int dim = mesh->Dimension();
185 
186  // 3. Refine the mesh to increase the resolution. In this example we do
187  // 'ref_levels' of uniform refinement. We choose 'ref_levels' to be the
188  // largest number that gives a final mesh with no more than 50,000
189  // elements.
190  {
191  if (ref_levels < 0)
192  {
193  ref_levels =
194  (int)floor(log(5000./mesh->GetNE())/log(2.)/dim);
195  }
196 
197  for (int l = 0; l < ref_levels; l++)
198  {
199  mesh->UniformRefinement();
200  }
201  mesh->PrintInfo();
202  }
203 
204  // 4. Define a finite element space on the mesh. Here we use continuous
205  // Lagrange finite elements of the specified order. If order < 1, we
206  // instead use an isoparametric/isogeometric space.
208  NURBSExtension *NURBSext = NULL;
209  int own_fec = 0;
210 
211  if (mesh->NURBSext)
212  {
213  fec = new NURBSFECollection(order[0]);
214  own_fec = 1;
215 
216  int nkv = mesh->NURBSext->GetNKV();
217  if (order.Size() == 1)
218  {
219  int tmp = order[0];
220  order.SetSize(nkv);
221  order = tmp;
222  }
223 
224  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
225  NURBSext = new NURBSExtension(mesh->NURBSext, order);
226 
227  // Read periodic BCs from file
228  std::ifstream in;
229  in.open(per_file, std::ifstream::in);
230  if (in.is_open())
231  {
232  int psize;
233  in >> psize;
234  master.SetSize(psize);
235  slave.SetSize(psize);
236  master.Load(in, psize);
237  slave.Load(in, psize);
238  in.close();
239  }
240  master.Print();
241  slave.Print();
242  NURBSext->ConnectBoundaries(master,slave);
243  }
244  else if (order[0] == -1) // Isoparametric
245  {
246  if (mesh->GetNodes())
247  {
248  fec = mesh->GetNodes()->OwnFEC();
249  own_fec = 0;
250  cout << "Using isoparametric FEs: " << fec->Name() << endl;
251  }
252  else
253  {
254  cout <<"Mesh does not have FEs --> Assume order 1.\n";
255  fec = new H1_FECollection(1, dim);
256  own_fec = 1;
257  }
258  }
259  else
260  {
261  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
262  fec = new H1_FECollection(abs(order[0]), dim);
263  own_fec = 1;
264  }
265 
266  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, NURBSext, fec);
267  cout << "Number of finite element unknowns: "
268  << fespace->GetTrueVSize() << endl;
269 
270  if (!ibp)
271  {
272  if (!mesh->NURBSext)
273  {
274  cout << "No integration by parts requires a NURBS mesh."<< endl;
275  return 2;
276  }
277  if (mesh->NURBSext->GetNP()>1)
278  {
279  cout << "No integration by parts requires a NURBS mesh, with only 1 patch."<<
280  endl;
281  cout << "A C_1 discretisation is required."<< endl;
282  cout << "Currently only C_0 multipatch coupling implemented."<< endl;
283  return 3;
284  }
285  if (order[0]<2)
286  {
287  cout << "No integration by parts requires at least quadratic NURBS."<< endl;
288  cout << "A C_1 discretisation is required."<< endl;
289  return 4;
290  }
291  }
292 
293  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
294  // In this example, the boundary conditions are defined by marking all
295  // the boundary attributes from the mesh as essential (Dirichlet) and
296  // converting them to a list of true dofs.
297  Array<int> ess_tdof_list;
298  if (mesh->bdr_attributes.Size())
299  {
300  Array<int> ess_bdr(mesh->bdr_attributes.Max());
301  if (strongBC)
302  {
303  ess_bdr = 1;
304  }
305  else
306  {
307  ess_bdr = 0;
308  }
309 
310  // Remove periodic BCs
311  for (int i = 0; i < master.Size(); i++)
312  {
313  ess_bdr[master[i]-1] = 0;
314  ess_bdr[slave[i]-1] = 0;
315  }
316  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
317  }
318 
319  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
320  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
321  // the basis functions in the finite element fespace.
322  ConstantCoefficient one(1.0);
323  ConstantCoefficient zero(0.0);
324 
325  LinearForm *b = new LinearForm(fespace);
327  if (!strongBC)
329  new DGDirichletLFIntegrator(zero, one, -1.0, kappa));
330  b->Assemble();
331 
332  // 7. Define the solution vector x as a finite element grid function
333  // corresponding to fespace. Initialize x with initial guess of zero,
334  // which satisfies the boundary conditions.
335  GridFunction x(fespace);
336  x = 0.0;
337 
338  // 8. Set up the bilinear form a(.,.) on the finite element space
339  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
340  // domain integrator.
341  BilinearForm *a = new BilinearForm(fespace);
342  if (ibp)
343  {
345  }
346  else
347  {
348  a->AddDomainIntegrator(new Diffusion2Integrator(one));
349  }
350 
351  if (!strongBC)
352  {
353  a->AddBdrFaceIntegrator(new DGDiffusionIntegrator(one, -1.0, kappa));
354  }
355 
356  // 9. Assemble the bilinear form and the corresponding linear system,
357  // applying any necessary transformations such as: eliminating boundary
358  // conditions, applying conforming constraints for non-conforming AMR,
359  // static condensation, etc.
360  if (static_cond) { a->EnableStaticCondensation(); }
361  a->Assemble();
362 
363  SparseMatrix A;
364  Vector B, X;
365  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
366 
367  cout << "Size of linear system: " << A.Height() << endl;
368 
369 #ifndef MFEM_USE_SUITESPARSE
370  // 10. Define a simple Jacobi preconditioner and use it to
371  // solve the system A X = B with PCG.
372  GSSmoother M(A);
373  PCG(A, M, B, X, 1, 200, 1e-12, 0.0);
374 #else
375  // 10. If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
376  UMFPackSolver umf_solver;
377  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
378  umf_solver.SetOperator(A);
379  umf_solver.Mult(B, X);
380 #endif
381 
382  // 11. Recover the solution as a finite element grid function.
383  a->RecoverFEMSolution(X, *b, x);
384 
385  // 12. Save the refined mesh and the solution. This output can be viewed later
386  // using GLVis: "glvis -m refined.mesh -g sol.gf".
387  ofstream mesh_ofs("refined.mesh");
388  mesh_ofs.precision(8);
389  mesh->Print(mesh_ofs);
390  ofstream sol_ofs("sol.gf");
391  sol_ofs.precision(8);
392  x.Save(sol_ofs);
393 
394  // 13. Send the solution by socket to a GLVis server.
395  if (visualization)
396  {
397  char vishost[] = "localhost";
398  int visport = 19916;
399  socketstream sol_sock(vishost, visport);
400  sol_sock.precision(8);
401  sol_sock << "solution\n" << *mesh << x << flush;
402  }
403 
404  // 14. Save data in the VisIt format
405  VisItDataCollection visit_dc("Example1", mesh);
406  visit_dc.RegisterField("solution", &x);
407  visit_dc.Save();
408 
409  // 15. Free the used memory.
410  delete a;
411  delete b;
412  delete fespace;
413  if (own_fec) { delete fec; }
414  delete mesh;
415 
416  return 0;
417 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
Arbitrary order non-uniform rational B-splines (NURBS) finite elements.
Definition: fe_coll.hpp:602
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:97
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
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:30
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
virtual void PrintInfo(std::ostream &os=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
Definition: mesh.hpp:1756
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_base.cpp:202
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:102
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
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:355
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_base.hpp:327
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)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:564
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:928
virtual void Save()
Save the collection and a VisIt root file.
int Space() const
Returns the type of FunctionSpace on the element.
Definition: fe_base.hpp:337
double kappa
Definition: ex24.cpp:54
Data type for Gauss-Seidel smoother of sparse matrix.
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:1025
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3619
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
constexpr char vishost[]
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:154
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9737
constexpr int visport
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
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
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:904
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:83
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:566
int Dimension() const
Definition: mesh.hpp:999
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
double Control[UMFPACK_CONTROL]
Definition: solvers.hpp:1036
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
int GetNP() const
Definition: nurbs.hpp:346
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
FDualNumber< tbase > log(const FDualNumber< tbase > &f)
log([dual number])
Definition: fdual.hpp:515
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:272
A &quot;square matrix&quot; operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
virtual const char * Name() const
Definition: fe_coll.hpp:61
Class for integration point with weight.
Definition: intrules.hpp:25
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:1646
void ConnectBoundaries()
Definition: nurbs.cpp:1793
int dim
Definition: ex24.cpp:53
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
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:60
IntegrationRules RefinedIntRules(1, Quadrature1D::GaussLegendre)
A global object with all refined integration rules.
Definition: intrules.hpp:382
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7841
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:23
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
int main()
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3179
void EnableStaticCondensation()
Enable the use of static condensation. For details see the description for class StaticCondensation i...
void AddBdrFaceIntegrator(BilinearFormIntegrator *bfi)
Adds new boundary Face Integrator. Assumes ownership of bfi.
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:3084
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150