MFEM  v4.1.0
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 square-nurbs.mesh -o 2 -no-ibp
26 // mpirun -np 4 nurbs_ex1p -m cube-nurbs.mesh -o 2 -no-ibp
27 // mpirun -np 4 nurbs_ex1p -m 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 j = 0; j < nd; j++)
125  {
126  for (int i = 0; i < nd; i++)
127  {
128  elmat(i, j) += w*shape(i)*laplace(j);
129  }
130  }
131  }
132  }
133 
134 };
135 
136 int main(int argc, char *argv[])
137 {
138  // 1. Initialize MPI.
139  int num_procs, myid;
140  MPI_Init(&argc, &argv);
141  MPI_Comm_size(MPI_COMM_WORLD, &num_procs);
142  MPI_Comm_rank(MPI_COMM_WORLD, &myid);
143 
144  // 2. Parse command-line options.
145  const char *mesh_file = "../../data/star.mesh";
146  Array<int> order(1);
147  order[0] = 1;
148  bool static_cond = false;
149  bool visualization = 1;
150  bool ibp = 1;
151 
152  OptionsParser args(argc, argv);
153  args.AddOption(&mesh_file, "-m", "--mesh",
154  "Mesh file to use.");
155  args.AddOption(&order, "-o", "--order",
156  "Finite element order (polynomial degree) or -1 for"
157  " isoparametric space.");
158  args.AddOption(&ibp, "-ibp", "--ibp", "-no-ibp",
159  "--no-ibp",
160  "Selects the standard weak form (IBP) or the nonstandard (NO-IBP).");
161  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
162  "--no-static-condensation", "Enable static condensation.");
163  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
164  "--no-visualization",
165  "Enable or disable GLVis visualization.");
166  args.Parse();
167  if (!args.Good())
168  {
169  if (myid == 0)
170  {
171  args.PrintUsage(cout);
172  }
173  MPI_Finalize();
174  return 1;
175  }
176  if (myid == 0)
177  {
178  args.PrintOptions(cout);
179  }
180 
181  // 3. Read the (serial) mesh from the given mesh file on all processors. We
182  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
183  // and volume meshes with the same code.
184  Mesh *mesh = new Mesh(mesh_file, 1, 1);
185  int dim = mesh->Dimension();
186 
187  // 4. Refine the serial mesh on all processors to increase the resolution. In
188  // this example we do 'ref_levels' of uniform refinement. We choose
189  // 'ref_levels' to be the largest number that gives a final mesh with no
190  // more than 10,000 elements.
191  {
192  int ref_levels =
193  (int)floor(log(10000./mesh->GetNE())/log(2.)/dim);
194  for (int l = 0; l < ref_levels; l++)
195  {
196  mesh->UniformRefinement();
197  }
198  }
199 
200  // 5. Define a parallel mesh by a partitioning of the serial mesh. Refine
201  // this mesh further in parallel to increase the resolution. Once the
202  // parallel mesh is defined, the serial mesh can be deleted.
203  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
204  delete mesh;
205  if (!pmesh->NURBSext)
206  {
207  int par_ref_levels = 2;
208  for (int l = 0; l < par_ref_levels; l++)
209  {
210  pmesh->UniformRefinement();
211  }
212  }
213 
214  // 6. Define a parallel finite element space on the parallel mesh. Here we
215  // use continuous Lagrange finite elements of the specified order. If
216  // order < 1, we instead use an isoparametric/isogeometric space.
218  NURBSExtension *NURBSext = NULL;
219  int own_fec = 0;
220 
221  if (order[0] == -1) // Isoparametric
222  {
223  if (pmesh->GetNodes())
224  {
225  fec = pmesh->GetNodes()->OwnFEC();
226  own_fec = 0;
227  cout << "Using isoparametric FEs: " << fec->Name() << endl;
228  }
229  else
230  {
231  cout <<"Mesh does not have FEs --> Assume order 1.\n";
232  fec = new H1_FECollection(1, dim);
233  own_fec = 1;
234  }
235  }
236  else if (pmesh->NURBSext && (order[0] > 0) ) // Subparametric NURBS
237  {
238  fec = new NURBSFECollection(order[0]);
239  own_fec = 1;
240  int nkv = pmesh->NURBSext->GetNKV();
241 
242  if (order.Size() == 1)
243  {
244  int tmp = order[0];
245  order.SetSize(nkv);
246  order = tmp;
247  }
248  if (order.Size() != nkv ) { mfem_error("Wrong number of orders set."); }
249  NURBSext = new NURBSExtension(pmesh->NURBSext, order);
250  }
251  else
252  {
253  if (order.Size() > 1) { cout <<"Wrong number of orders set, needs one.\n"; }
254  fec = new H1_FECollection(abs(order[0]), dim);
255  own_fec = 1;
256  }
257  ParFiniteElementSpace *fespace = new ParFiniteElementSpace(pmesh,NURBSext,fec);
258  HYPRE_Int size = fespace->GlobalTrueVSize();
259  if (myid == 0)
260  {
261  cout << "Number of finite element unknowns: " << size << endl;
262  }
263 
264  if (!ibp)
265  {
266  if (!pmesh->NURBSext)
267  {
268  cout << "No integration by parts requires a NURBS mesh."<< endl;
269  return 2;
270  }
271  if (pmesh->NURBSext->GetNP()>1)
272  {
273  cout << "No integration by parts requires a NURBS mesh, with only 1 patch."<<
274  endl;
275  cout << "A C_1 discretisation is required."<< endl;
276  cout << "Currently only C_0 multipatch coupling implemented."<< endl;
277  return 3;
278  }
279  if (order[0]<2)
280  {
281  cout << "No integration by parts requires at least quadratic NURBS."<< endl;
282  cout << "A C_1 discretisation is required."<< endl;
283  return 4;
284  }
285  }
286 
287  // 7. Determine the list of true (i.e. parallel conforming) essential
288  // boundary dofs. In this example, the boundary conditions are defined
289  // by marking all the boundary attributes from the mesh as essential
290  // (Dirichlet) and converting them to a list of true dofs.
291  Array<int> ess_tdof_list;
292  if (pmesh->bdr_attributes.Size())
293  {
294  Array<int> ess_bdr(pmesh->bdr_attributes.Max());
295  ess_bdr = 1;
296  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
297  }
298 
299  // 8. Set up the parallel linear form b(.) which corresponds to the
300  // right-hand side of the FEM linear system, which in this case is
301  // (1,phi_i) where phi_i are the basis functions in fespace.
302  ParLinearForm *b = new ParLinearForm(fespace);
303  ConstantCoefficient one(1.0);
305  b->Assemble();
306 
307  // 9. Define the solution vector x as a parallel finite element grid function
308  // corresponding to fespace. Initialize x with initial guess of zero,
309  // which satisfies the boundary conditions.
310  ParGridFunction x(fespace);
311  x = 0.0;
312 
313  // 10. Set up the parallel bilinear form a(.,.) on the finite element space
314  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
315  // domain integrator.
316  ParBilinearForm *a = new ParBilinearForm(fespace);
317  if (ibp)
318  {
320  }
321  else
322  {
323  a->AddDomainIntegrator(new Diffusion2Integrator(one));
324  }
325 
326  // 11. Assemble the parallel bilinear form and the corresponding linear
327  // system, applying any necessary transformations such as: parallel
328  // assembly, eliminating boundary conditions, applying conforming
329  // constraints for non-conforming AMR, static condensation, etc.
330  if (static_cond) { a->EnableStaticCondensation(); }
331  a->Assemble();
332 
333  HypreParMatrix A;
334  Vector B, X;
335  a->FormLinearSystem(ess_tdof_list, x, *b, A, X, B);
336 
337  if (myid == 0)
338  {
339  cout << "Size of linear system: " << A.GetGlobalNumRows() << endl;
340  }
341 
342  // 12. Define and apply a parallel PCG solver for AX=B with the BoomerAMG
343  // preconditioner from hypre.
344  HypreSolver *amg = new HypreBoomerAMG(A);
345  HyprePCG *pcg = new HyprePCG(A);
346  pcg->SetTol(1e-12);
347  pcg->SetMaxIter(200);
348  pcg->SetPrintLevel(2);
349  pcg->SetPreconditioner(*amg);
350  pcg->Mult(B, X);
351 
352  // 13. Recover the parallel grid function corresponding to X. This is the
353  // local finite element solution on each processor.
354  a->RecoverFEMSolution(X, *b, x);
355 
356  // 14. Save the refined mesh and the solution in parallel. This output can
357  // be viewed later using GLVis: "glvis -np <np> -m mesh -g sol".
358  {
359  ostringstream mesh_name, sol_name;
360  mesh_name << "mesh." << setfill('0') << setw(6) << myid;
361  sol_name << "sol." << setfill('0') << setw(6) << myid;
362 
363  ofstream mesh_ofs(mesh_name.str().c_str());
364  mesh_ofs.precision(8);
365  pmesh->Print(mesh_ofs);
366 
367  ofstream sol_ofs(sol_name.str().c_str());
368  sol_ofs.precision(8);
369  x.Save(sol_ofs);
370  }
371 
372  // 15. Send the solution by socket to a GLVis server.
373  if (visualization)
374  {
375  char vishost[] = "localhost";
376  int visport = 19916;
377  socketstream sol_sock(vishost, visport);
378  sol_sock << "parallel " << num_procs << " " << myid << "\n";
379  sol_sock.precision(8);
380  sol_sock << "solution\n" << *pmesh << x << flush;
381  }
382 
383  // 16. Save data in the VisIt format
384  VisItDataCollection visit_dc("Example1-Parallel", pmesh);
385  visit_dc.RegisterField("solution", &x);
386  visit_dc.Save();
387 
388  // 17. Free the used memory.
389  delete pcg;
390  delete amg;
391  delete a;
392  delete b;
393  delete fespace;
394  if (own_fec) { delete fec; }
395  delete pmesh;
396 
397  MPI_Finalize();
398 
399  return 0;
400 }
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
void SetTol(double tol)
Definition: hypre.cpp:2305
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 GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:770
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:308
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:890
Subclass constant coefficient.
Definition: coefficient.hpp:67
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(...
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
int GetNKV() const
Definition: nurbs.hpp:350
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
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:693
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:485
int Space() const
Returns the type of space on each element.
Definition: fe.hpp:328
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:62
void SetPrintLevel(int print_lvl)
Definition: hypre.cpp:2320
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe.hpp:311
The BoomerAMG solver in hypre.
Definition: hypre.hpp:951
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_Int GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:421
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
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4102
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetMaxIter(int max_iter)
Definition: hypre.cpp:2310
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
PCG solver in hypre.
Definition: hypre.hpp:743
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 RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
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
int dim
Definition: ex24.cpp:43
void SetPreconditioner(HypreSolver &precond)
Set the hypre solver to be used as a preconditioner.
Definition: hypre.cpp:2325
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
Class for parallel bilinear form.
Abstract class for hypre&#39;s solvers and preconditioners.
Definition: hypre.hpp:684
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 parallel grid function.
Definition: pgridfunc.hpp:32
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:90
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:187
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve Ax=b with hypre&#39;s PCG.
Definition: hypre.cpp:2348
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:376
void EnableStaticCondensation()
bool Good() const
Definition: optparser.hpp:125