MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
nurbs_ex1p.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - Parallel NURBS Version
2 //
3 // Compile with: make nurbs_ex1p
4 //
5 // Sample runs: mpirun -np 4 nurbs_ex1p -m ../../data/square-disc.mesh
6 // mpirun -np 4 nurbs_ex1p -m ../../data/star.mesh
7 // mpirun -np 4 nurbs_ex1p -m ../../data/escher.mesh
8 // mpirun -np 4 nurbs_ex1p -m ../../data/fichera.mesh
9 // mpirun -np 4 nurbs_ex1p -m ../../data/square-disc-p2.vtk -o 2
10 // mpirun -np 4 nurbs_ex1p -m ../../data/square-disc-p3.mesh -o 3
11 // mpirun -np 4 nurbs_ex1p -m ../../data/square-disc-nurbs.mesh -o -1
12 // mpirun -np 4 nurbs_ex1p -m ../../data/disc-nurbs.mesh -o -1
13 // mpirun -np 4 nurbs_ex1p -m ../../data/pipe-nurbs.mesh -o -1
14 // mpirun -np 4 nurbs_ex1p -m ../../data/ball-nurbs.mesh -o 2
15 // mpirun -np 4 nurbs_ex1p -m ../../data/star-surf.mesh
16 // mpirun -np 4 nurbs_ex1p -m ../../data/square-disc-surf.mesh
17 // mpirun -np 4 nurbs_ex1p -m ../../data/inline-segment.mesh
18 // mpirun -np 4 nurbs_ex1p -m ../../data/amr-quad.mesh
19 // mpirun -np 4 nurbs_ex1p -m ../../data/amr-hex.mesh
20 // mpirun -np 4 nurbs_ex1p -m ../../data/mobius-strip.mesh
21 // mpirun -np 4 nurbs_ex1p -m ../../data/mobius-strip.mesh -o -1 -sc
22 // mpirun -np 4 nurbs_ex1p -m ../../data/square-disc-nurbs.mesh -o -1
23 // mpirun -np 4 nurbs_ex1p -m ../../data/disc-nurbs.mesh -o -1
24 // mpirun -np 4 nurbs_ex1p -m ../../data/pipe-nurbs.mesh -o -1
25 // mpirun -np 4 nurbs_ex1p -m ../../data/square-nurbs.mesh -o 2 -no-ibp
26 // mpirun -np 4 nurbs_ex1p -m ../../data/cube-nurbs.mesh -o 2 -no-ibp
27 // mpirun -np 4 nurbs_ex1p -m ../../data/pipe-nurbs-2d.mesh -o 2 -no-ibp
28 
29 // Description: This example code demonstrates the use of MFEM to define a
30 // simple finite element discretization of the Laplace problem
31 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
32 // Specifically, we discretize using a FE space of the specified
33 // order, or if order < 1 using an isoparametric/isogeometric
34 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
35 // NURBS mesh, etc.)
36 //
37 // The example highlights the use of mesh refinement, finite
38 // element grid functions, as well as linear and bilinear forms
39 // corresponding to the left-hand side and right-hand side of the
40 // discrete linear system. We also cover the explicit elimination
41 // of essential boundary conditions, static condensation, and the
42 // optional connection to the GLVis tool for visualization.
43 
44 #include "mfem.hpp"
45 #include <fstream>
46 #include <iostream>
47 
48 using namespace std;
49 using namespace mfem;
50 /** Class for integrating the bilinear form a(u,v) := (Q Laplace u, v) where Q
51  can be a scalar coefficient. */
52 class Diffusion2Integrator: public BilinearFormIntegrator
53 {
54 private:
55 #ifndef MFEM_THREAD_SAFE
56  Vector shape,laplace;
57 #endif
58  Coefficient *Q;
59 
60 public:
61  /// Construct a diffusion integrator with coefficient Q = 1
62  Diffusion2Integrator() { Q = NULL; }
63 
64  /// Construct a diffusion integrator with a scalar coefficient q
65  Diffusion2Integrator (Coefficient &q) : Q(&q) { }
66 
67  /** Given a particular Finite Element
68  computes the element stiffness matrix elmat. */
69  virtual void AssembleElementMatrix(const FiniteElement &el,
71  DenseMatrix &elmat)
72  {
73  int nd = el.GetDof();
74  int dim = el.GetDim();
75  double w;
76 
77 #ifdef MFEM_THREAD_SAFE
78  Vector shape(nd);
79  Vector laplace(nd);
80 #else
81  shape.SetSize(nd);
82  laplace.SetSize(nd);
83 #endif
84  elmat.SetSize(nd);
85 
86  const IntegrationRule *ir = IntRule;
87  if (ir == NULL)
88  {
89  int order;
90  if (el.Space() == FunctionSpace::Pk)
91  {
92  order = 2*el.GetOrder() - 2;
93  }
94  else
95  {
96  order = 2*el.GetOrder() + dim - 1;
97  }
98 
99  if (el.Space() == FunctionSpace::rQk)
100  {
101  ir = &RefinedIntRules.Get(el.GetGeomType(),order);
102  }
103  else
104  {
105  ir = &IntRules.Get(el.GetGeomType(),order);
106  }
107  }
108 
109  elmat = 0.0;
110  for (int i = 0; i < ir->GetNPoints(); i++)
111  {
112  const IntegrationPoint &ip = ir->IntPoint(i);
113  Trans.SetIntPoint(&ip);
114  w = -ip.weight * Trans.Weight();
115 
116  el.CalcShape(ip, shape);
117  el.CalcPhysLaplacian(Trans, laplace);
118 
119  if (Q)
120  {
121  w *= Q->Eval(Trans, ip);
122  }
123 
124  for (int jj = 0; jj < nd; jj++)
125  {
126  for (int ii = 0; ii < nd; ii++)
127  {
128  elmat(ii, jj) += w*shape(ii)*laplace(jj);
129  }
130  }
131  }
132  }
133 
134 };
135 
136 int main(int argc, char *argv[])
137 {
138  // 1. Initialize MPI and HYPRE.
139  Mpi::Init(argc, argv);
140  int num_procs = Mpi::WorldSize();
141  int myid = Mpi::WorldRank();
142  Hypre::Init();
143 
144  // 2. Parse command-line options.
145  const char *mesh_file = "../../data/star.mesh";
146  int ref_levels = -1;
147  Array<int> order(1);
148  order[0] = 1;
149  bool static_cond = false;
150  bool visualization = 1;
151  bool ibp = 1;
152  bool strongBC = 1;
153  double kappa = -1;
154 
155  OptionsParser args(argc, argv);
156  args.AddOption(&mesh_file, "-m", "--mesh",
157  "Mesh file to use.");
158  args.AddOption(&ref_levels, "-r", "--refine",
159  "Number of times to refine the mesh uniformly, -1 for auto.");
160  args.AddOption(&order, "-o", "--order",
161  "Finite element order (polynomial degree) or -1 for"
162  " isoparametric space.");
163  args.AddOption(&ibp, "-ibp", "--ibp", "-no-ibp",
164  "--no-ibp",
165  "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
166  args.AddOption(&strongBC, "-sbc", "--strong-bc", "-wbc",
167  "--weak-bc",
168  "Selects strong or weak enforcement of Dirichlet BCs.");
169  args.AddOption(&kappa, "-k", "--kappa",
170  "Sets the SIPG penalty parameters, should be positive."
171  " Negative values are replaced with (order+1)^2.");
172  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
173  "--no-static-condensation", "Enable static condensation.");
174  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
175  "--no-visualization",
176  "Enable or disable GLVis visualization.");
177  args.Parse();
178  if (!args.Good())
179  {
180  if (myid == 0)
181  {
182  args.PrintUsage(cout);
183  }
184  return 1;
185  }
186  if (!strongBC & (kappa < 0))
187  {
188  kappa = 4*(order.Max()+1)*(order.Max()+1);
189  }
190  if (myid == 0)
191  {
192  args.PrintOptions(cout);
193  }
194 
195  // 3. Read the (serial) mesh from the given mesh file on all processors. We
196  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
197  // and volume meshes with the same code.
198  Mesh *mesh = new Mesh(mesh_file, 1, 1);
199  int dim = mesh->Dimension();
200 
201  // 4. Refine the serial mesh on all processors to increase the resolution. In
202  // this example we do 'ref_levels' of uniform refinement. We choose
203  // 'ref_levels' to be the largest number that gives a final mesh with no
204  // more than 10,000 elements.
205  {
206  if (ref_levels < 0)
207  {
208  ref_levels =
209  (int)floor(log(5000./mesh->GetNE())/log(2.)/dim);
210  }
211 
212  for (int l = 0; l < ref_levels; l++)
213  {
214  mesh->UniformRefinement();
215  }
216  if (myid == 0)
217  {
218  mesh->PrintInfo();
219  }
220  }
221 
222  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
223  // this mesh further in parallel to increase the resolution. Once the
224  // parallel mesh is defined, the serial mesh can be deleted.
225  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
226  delete mesh;
227  if (!pmesh->NURBSext)
228  {
229  int par_ref_levels = 2;
230  for (int l = 0; l < par_ref_levels; l++)
231  {
232  pmesh->UniformRefinement();
233  }
234  }
235 
236  // 6. Define a parallel finite element space on the parallel mesh. Here we
237  // use continuous Lagrange finite elements of the specified order. If
238  // order < 1, we instead use an isoparametric/isogeometric space.
240  NURBSExtension *NURBSext = NULL;
241  int own_fec = 0;
242 
243  if (order[0] == -1) // Isoparametric
244  {
245  if (pmesh->GetNodes())
246  {
247  fec = pmesh->GetNodes()->OwnFEC();
248  own_fec = 0;
249  cout << "Using isoparametric FEs: " << fec->Name() << endl;
250  }
251  else
252  {
253  cout <<"Mesh does not have FEs --> Assume order 1.\n";
254  fec = new H1_FECollection(1, dim);
255  own_fec = 1;
256  }
257  }
258  else if (pmesh->NURBSext && (order[0] > 0) ) // Subparametric NURBS
259  {
260  fec = new NURBSFECollection(order[0]);
261  own_fec = 1;
262  int nkv = pmesh->NURBSext->GetNKV();
263 
264  if (order.Size() == 1)
265  {
266  int tmp = order[0];
267  order.SetSize(nkv);
268  order = tmp;
269  }
270  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
271  NURBSext = new NURBSExtension(pmesh->NURBSext, order);
272  }
273  else
274  {
275  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
276  fec = new H1_FECollection(abs(order[0]), dim);
277  own_fec = 1;
278  }
279  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh,NURBSext,fec);
280  HYPRE_BigInt size = fespace->GlobalTrueVSize();
281  if (myid == 0)
282  {
283  cout << "Number of finite element unknowns: " << size << endl;
284  }
285 
286  if (!ibp)
287  {
288  if (!pmesh->NURBSext)
289  {
290  cout << "No integration by parts requires a NURBS mesh."<< endl;
291  return 2;
292  }
293  if (pmesh->NURBSext->GetNP()>1)
294  {
295  cout << "No integration by parts requires a NURBS mesh, with only 1 patch."<<
296  endl;
297  cout << "A C_1 discretisation is required."<< endl;
298  cout << "Currently only C_0 multipatch coupling implemented."<< endl;
299  return 3;
300  }
301  if (order[0]<2)
302  {
303  cout << "No integration by parts requires at least quadratic NURBS."<< endl;
304  cout << "A C_1 discretisation is required."<< endl;
305  return 4;
306  }
307  }
308 
309  // 7. Determine the list of true (i.e. parallel conforming) essential
310  // boundary dofs. In this example, the boundary conditions are defined
311  // by marking all the boundary attributes from the mesh as essential
312  // (Dirichlet) and converting them to a list of true dofs.
313  Array<int> ess_tdof_list;
314  if (pmesh->bdr_attributes.Size())
315  {
316  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
317  if (strongBC)
318  {
319  ess_bdr = 1;
320  }
321  else
322  {
323  ess_bdr = 0;
324  }
325  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
326  }
327 
328  // 8. Set up the parallel linear form b(.) which corresponds to the
329  // right-hand side of the FEM linear system, which in this case is
330  // (1,phi_i) where phi_i are the basis functions in fespace.
331  ConstantCoefficient one(1.0);
332  ConstantCoefficient zero(0.0);
333 
334  ParLinearForm *b = new ParLinearForm(fespace);
336 
337  if (!strongBC)
339  new DGDirichletLFIntegrator(zero, one, -1.0, kappa));
340  b->Assemble();
341 
342  // 9. Define the solution vector x as a parallel finite element grid function
343  // corresponding to fespace. Initialize x with initial guess of zero,
344  // which satisfies the boundary conditions.
345  ParGridFunction x(fespace);
346  x = 0.0;
347 
348  // 10. Set up the parallel bilinear form a(.,.) on the finite element space
349  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
350  // domain integrator.
351  ParBilinearForm *a = new ParBilinearForm(fespace);
352  if (ibp)
353  {
355  }
356  else
357  {
358  a->AddDomainIntegrator(new Diffusion2Integrator(one));
359  }
360  if (!strongBC)
361  {
362  a->AddBdrFaceIntegrator(new DGDiffusionIntegrator(one, -1.0, kappa));
363  }
364 
365  // 11. Assemble the parallel bilinear form and the corresponding linear
366  // system, applying any necessary transformations such as: parallel
367  // assembly, eliminating boundary conditions, applying conforming
368  // constraints for non-conforming AMR, static condensation, etc.
369  if (static_cond) { a->EnableStaticCondensation(); }
370  a->Assemble();
371 
372  HypreParMatrix A;
373  Vector B, X;
374  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
375 
376  if (myid == 0)
377  {
378  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
379  }
380 
381  // 12. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
382  // preconditioner from hypre.
383  HypreSolver *amg = new HypreBoomerAMG(A);
384  HyprePCG *pcg = new HyprePCG(A);
385  pcg->SetTol(1e-12);
386  pcg->SetMaxIter(200);
387  pcg->SetPrintLevel(2);
388  pcg->SetPreconditioner(*amg);
389  pcg->Mult(B, X);
390 
391  // 13. Recover the parallel grid function corresponding to X. This is the
392  // local finite element solution on each processor.
393  a->RecoverFEMSolution(X, *b, x);
394 
395  // 14. Save the refined mesh and the solution in parallel. This output can
396  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
397  {
398  ostringstream mesh_name, sol_name;
399  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
400  sol_name << "sol." << setfill('0') << setw(6) << myid;
401 
402  ofstream mesh_ofs(mesh_name.str().c_str());
403  mesh_ofs.precision(8);
404  pmesh->Print(mesh_ofs);
405 
406  ofstream sol_ofs(sol_name.str().c_str());
407  sol_ofs.precision(8);
408  x.Save(sol_ofs);
409  }
410 
411  // 15. Send the solution by socket to a GLVis server.
412  if (visualization)
413  {
414  char vishost[] = "localhost";
415  int visport = 19916;
416  socketstream sol_sock(vishost, visport);
417  sol_sock << "parallel " << num_procs << " " << myid << "\n";
418  sol_sock.precision(8);
419  sol_sock << "solution\n" << *pmesh << x << flush;
420  }
421 
422  // 16. Save data in the VisIt format
423  VisItDataCollection visit_dc("Example1-Parallel", pmesh);
424  visit_dc.RegisterField("solution", &x);
425  visit_dc.Save();
426 
427  // 17. Free the used memory.
428  delete pcg;
429  delete amg;
430  delete a;
431  delete b;
432  delete fespace;
433  if (own_fec) { delete fec; }
434  delete pmesh;
435 
436  return 0;
437 }
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:606
void SetTol(double tol)
Definition: hypre.cpp:3995
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:108
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
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(...
virtual void PrintInfo(std::ostream &os=mfem::out)
In serial, this method calls PrintCharacteristics(). In parallel, additional information about the pa...
Definition: mesh.hpp:1781
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 SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
HYPRE_BigInt GlobalTrueVSize() const
Definition: pfespace.hpp:285
int GetNKV() const
Definition: nurbs.hpp:355
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
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:873
int Space() const
Returns the type of FunctionSpace on the element.
Definition: fe_base.hpp:337
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:4015
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1579
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
Class for parallel linear form.
Definition: plinearform.hpp:26
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
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
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Definition: plinearform.cpp:46
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
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
void Assemble(int skip_zeros=1)
Assemble the local matrix.
double Weight()
Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:131
void AddBdrFaceIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Face Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:85
int Dimension() const
Definition: mesh.hpp:1006
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
void SetMaxIter(int max_iter)
Definition: hypre.cpp:4005
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:41
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:35
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
PCG solver in hypre.
Definition: hypre.hpp:1204
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: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
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
HYPRE_Int HYPRE_BigInt
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:272
virtual const char * Name() const
Definition: fe_coll.hpp:65
Class for integration point with weight.
Definition: intrules.hpp:25
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:635
int dim
Definition: ex24.cpp:53
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:4020
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
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:1092
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:7908
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
void Print(std::ostream &out=mfem::out) const override
Definition: pmesh.cpp:4770
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:4043
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
int main()
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 RegisterField(const std::string &field_name, GridFunction *gf) override
Add a grid function to the collection and update the root file.
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150