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