MFEM  v4.6.0
Finite element discretization library
nurbs_patch_ex1.cpp
Go to the documentation of this file.
1 // MFEM Example 1 - NURBS with patch-wise assembly
2 //
3 // Compile with: make nurbs_patch_ex1
4 //
5 // Sample runs: nurbs_patch_ex1 -incdeg 3 -ref 2 -iro 8 -patcha
6 // nurbs_patch_ex1 -incdeg 3 -ref 2 -iro 8 -patcha -pa
7 // nurbs_patch_ex1 -incdeg 3 -ref 2 -iro 8 -patcha -fint
8 //
9 // Description: This example code demonstrates the use of MFEM to define a
10 // simple finite element discretization of the Laplace problem
11 // -Delta u = 1 with homogeneous Dirichlet boundary conditions.
12 // Specifically, we discretize using a FE space of the specified
13 // order, or if order < 1 using an isoparametric/isogeometric
14 // space (i.e. quadratic for quadratic curvilinear mesh, NURBS for
15 // NURBS mesh, etc.)
16 //
17 // This example is a specialization of ex1 which demonstrates
18 // patch-wise matrix assembly and partial assembly on NURBS
19 // meshes. There is the option to compare run times of patch
20 // and element assembly, as well as relative error computation.
21 
22 #include "mfem.hpp"
23 #include <fstream>
24 #include <iostream>
25 
26 using namespace std;
27 using namespace mfem;
28 
30  Array<int> const& ess_tdof_list, const bool pa,
31  const bool algebraic_ceed, GridFunction & x);
32 
33 int main(int argc, char *argv[])
34 {
35  // 1. Parse command-line options.
36  const char *mesh_file = "../../data/beam-hex-nurbs.mesh";
37  int order = -1;
38  bool pa = false;
39  const char *device_config = "cpu";
40  bool visualization = true;
41  bool algebraic_ceed = false;
42  bool patchAssembly = false;
43  bool reducedIntegration = true;
44  bool compareToElementWise = true;
45  int nurbs_degree_increase = 0; // Elevate the NURBS mesh degree by this
46  int ref_levels = 0;
47  int ir_order = -1;
48 
49  OptionsParser args(argc, argv);
50  args.AddOption(&mesh_file, "-m", "--mesh",
51  "Mesh file to use.");
52  args.AddOption(&pa, "-pa", "--partial-assembly", "-no-pa",
53  "--no-partial-assembly", "Enable Partial Assembly.");
54  args.AddOption(&device_config, "-d", "--device",
55  "Device configuration string, see Device::Configure().");
56 #ifdef MFEM_USE_CEED
57  args.AddOption(&algebraic_ceed, "-a", "--algebraic", "-no-a", "--no-algebraic",
58  "Use algebraic Ceed solver");
59 #endif
60  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
61  "--no-visualization",
62  "Enable or disable GLVis visualization.");
63  args.AddOption(&patchAssembly, "-patcha", "--patch-assembly", "-no-patcha",
64  "--no-patch-assembly", "Enable patch-wise assembly.");
65  args.AddOption(&reducedIntegration, "-rint", "--reduced-integration", "-fint",
66  "--full-integration", "Enable reduced integration rules.");
67  args.AddOption(&ref_levels, "-ref", "--refine",
68  "Number of uniform mesh refinements.");
69  args.AddOption(&ir_order, "-iro", "--integration-order",
70  "Order of integration rule.");
71  args.AddOption(&nurbs_degree_increase, "-incdeg", "--nurbs-degree-increase",
72  "Elevate NURBS mesh degree by this amount.");
73  args.AddOption(&compareToElementWise, "-cew", "--compare-element",
74  "-no-compare", "-no-compare-element",
75  "Compute element-wise solution for comparison");
76  args.Parse();
77  if (!args.Good())
78  {
79  args.PrintUsage(cout);
80  return 1;
81  }
82  args.PrintOptions(cout);
83 
84  MFEM_VERIFY(!(pa && !patchAssembly), "Patch assembly must be used with -pa");
85 
86  // 2. Enable hardware devices such as GPUs, and programming models such as
87  // CUDA, OCCA, RAJA and OpenMP based on command line options.
88  Device device(device_config);
89  device.Print();
90 
91  // 3. Read the mesh from the given mesh file. For this NURBS patch example,
92  // only 3D hexahedral meshes are currently supported. The NURBS degree is
93  // optionally increased.
94  Mesh mesh(mesh_file, 1, 1);
95  int dim = mesh.Dimension();
96 
97  if (nurbs_degree_increase > 0) { mesh.DegreeElevate(nurbs_degree_increase); }
98 
99  // 4. Refine the mesh to increase the resolution.
100  for (int l = 0; l < ref_levels; l++)
101  {
102  mesh.UniformRefinement();
103  }
104 
105  // 5. Define an isoparametric/isogeometric finite element space on the mesh.
106  FiniteElementCollection *fec = nullptr;
107  bool delete_fec;
108  if (mesh.GetNodes())
109  {
110  fec = mesh.GetNodes()->OwnFEC();
111  delete_fec = false;
112  cout << "Using isoparametric FEs: " << fec->Name() << endl;
113  }
114  else
115  {
116  MFEM_ABORT("Mesh must have nodes");
117  }
118  FiniteElementSpace fespace(&mesh, fec);
119  cout << "Number of finite element unknowns: "
120  << fespace.GetTrueVSize() << endl;
121 
122  // 6. Determine the list of true (i.e. conforming) essential boundary dofs.
123  // In this example, the boundary conditions are defined by marking all
124  // the boundary attributes from the mesh as essential (Dirichlet) and
125  // converting them to a list of true dofs.
126  Array<int> ess_tdof_list;
127  if (mesh.bdr_attributes.Size())
128  {
129  Array<int> ess_bdr(mesh.bdr_attributes.Max());
130  ess_bdr = 1;
131  fespace.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
132  }
133 
134  // 7. Set up the linear form b(.) which corresponds to the right-hand side of
135  // the FEM linear system, which in this case is (1,phi_i) where phi_i are
136  // the basis functions in the finite element fespace.
137  LinearForm b(&fespace);
138  ConstantCoefficient one(1.0);
139  b.AddDomainIntegrator(new DomainLFIntegrator(one));
140  b.Assemble();
141 
142  // 8. Define the solution vector x as a finite element grid function
143  // corresponding to fespace. Initialize x with initial guess of zero,
144  // which satisfies the boundary conditions.
145  GridFunction x(&fespace);
146  x = 0.0;
147 
148  // 9. Set up the bilinear form a(.,.) on the finite element space
149  // corresponding to the Laplacian operator -Delta, by adding the Diffusion
150  // domain integrator.
152 
153  if (patchAssembly && reducedIntegration && !pa)
154  {
155  di->SetIntegrationMode(NonlinearFormIntegrator::Mode::PATCHWISE_REDUCED);
156  }
157  else if (patchAssembly)
158  {
159  di->SetIntegrationMode(NonlinearFormIntegrator::Mode::PATCHWISE);
160  }
161 
162  NURBSMeshRules *patchRule = nullptr;
163  if (order < 0)
164  {
165  if (ir_order == -1) { ir_order = 2*fec->GetOrder(); }
166  cout << "Using ir_order " << ir_order << endl;
167 
168  patchRule = new NURBSMeshRules(mesh.NURBSext->GetNP(), dim);
169  // Loop over patches and set a different rule for each patch.
170  for (int p=0; p<mesh.NURBSext->GetNP(); ++p)
171  {
173  mesh.NURBSext->GetPatchKnotVectors(p, kv);
174 
175  std::vector<const IntegrationRule*> ir1D(dim);
176  const IntegrationRule *ir = &IntRules.Get(Geometry::SEGMENT, ir_order);
177 
178  // Construct 1D integration rules by applying the rule ir to each
179  // knot span.
180  for (int i=0; i<dim; ++i)
181  {
182  ir1D[i] = ir->ApplyToKnotIntervals(*kv[i]);
183  }
184 
185  patchRule->SetPatchRules1D(p, ir1D);
186  } // loop (p) over patches
187 
188  patchRule->Finalize(mesh);
189  di->SetNURBSPatchIntRule(patchRule);
190  }
191 
192  // 10. Assemble and solve the linear system
193  cout << "Assembling system patch-wise and solving" << endl;
194  AssembleAndSolve(b, di, ess_tdof_list, pa, algebraic_ceed, x);
195 
196  delete patchRule;
197 
198  // 11. Save the refined mesh and the solution. This output can be viewed
199  // later using GLVis: "glvis -m refined.mesh -g sol.gf".
200  ofstream mesh_ofs("refined.mesh");
201  mesh_ofs.precision(8);
202  mesh.Print(mesh_ofs);
203  ofstream sol_ofs("sol.gf");
204  sol_ofs.precision(8);
205  x.Save(sol_ofs);
206 
207  // 12. Send the solution by socket to a GLVis server.
208  if (visualization)
209  {
210  char vishost[] = "localhost";
211  int visport = 19916;
212  socketstream sol_sock(vishost, visport);
213  sol_sock.precision(8);
214  sol_sock << "solution\n" << mesh << x << flush;
215  }
216 
217  // 13. Optionally assemble element-wise and solve the linear system, to
218  // compare timings and compute relative error.
219  if (compareToElementWise)
220  {
221  Vector x_pw, x_ew;
222  x.GetTrueDofs(x_pw);
223 
224  cout << "Assembling system element-wise and solving" << endl;
226  // Element-wise partial assembly is not supported on NURBS meshes, so we
227  // pass pa = false here.
228  AssembleAndSolve(b, d, ess_tdof_list, false, algebraic_ceed, x);
229 
230  x.GetTrueDofs(x_ew);
231 
232  const double solNorm = x_ew.Norml2();
233  x_ew -= x_pw;
234 
235  cout << "Element-wise solution norm " << solNorm << endl;
236  cout << "Relative error of patch-wise solution "
237  << x_ew.Norml2() / solNorm << endl;
238  }
239 
240  // 14. Free the used memory.
241  if (delete_fec)
242  {
243  delete fec;
244  }
245 
246  return 0;
247 }
248 
249 // This function deletes bfi when the BilinearForm goes out of scope.
251  Array<int> const& ess_tdof_list, const bool pa,
252  const bool algebraic_ceed, GridFunction & x)
253 {
254  FiniteElementSpace *fespace = b.FESpace();
255  BilinearForm a(fespace);
256  if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
257 
258  a.AddDomainIntegrator(bfi); // Takes ownership of bfi
259 
260  StopWatch sw;
261  sw.Start();
262 
263  // Assemble the bilinear form and the corresponding linear system, applying
264  // any necessary transformations such as: eliminating boundary conditions,
265  // applying conforming constraints for non-conforming AMR, etc.
266  a.Assemble();
267 
268  sw.Stop();
269 
270  const double timeAssemble = sw.RealTime();
271 
272  sw.Clear();
273  sw.Start();
274 
275  OperatorPtr A;
276  Vector B, X;
277  a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
278 
279  sw.Stop();
280 
281  const double timeFormLinearSystem = sw.RealTime();
282 
283  cout << "Timing for Assemble: " << timeAssemble << " seconds" << endl;
284  cout << "Timing for FormLinearSystem: " << timeFormLinearSystem << " seconds"
285  << endl;
286  cout << "Timing for entire setup: " << timeAssemble + timeFormLinearSystem
287  << " seconds" << endl;
288 
289  sw.Clear();
290  sw.Start();
291 
292  // Solve the linear system A X = B.
293  if (!pa)
294  {
295 #ifndef MFEM_USE_SUITESPARSE
296  // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
297  GSSmoother M((SparseMatrix&)(*A));
298  PCG(*A, M, B, X, 1, 200, 1e-20, 0.0);
299 #else
300  // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
301  UMFPackSolver umf_solver;
302  umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
303  umf_solver.SetOperator(*A);
304  umf_solver.Mult(B, X);
305 #endif
306  }
307  else
308  {
309  if (UsesTensorBasis(*fespace))
310  {
311  if (algebraic_ceed)
312  {
313  ceed::AlgebraicSolver M(a, ess_tdof_list);
314  PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
315  }
316  else
317  {
318  OperatorJacobiSmoother M(a, ess_tdof_list);
319  PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
320  }
321  }
322  else
323  {
324  CG(*A, B, X, 1, 400, 1e-20, 0.0);
325  }
326  }
327 
328  sw.Stop();
329  cout << "Timing for solve " << sw.RealTime() << endl;
330 
331  // Recover the solution as a finite element grid function.
332  a.RecoverFEMSolution(X, b, x);
333 }
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
int visport
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
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition: fe_coll.hpp:225
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
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
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition: device.cpp:279
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition: tic_toc.cpp:429
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
void Finalize(Mesh const &mesh)
Finalize() must be called before this class can be used for assembly. In particular, it defines data used by GetPointElement().
Definition: intrules.cpp:1925
STL namespace.
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition: fespace.hpp:1306
int main(int argc, char *argv[])
IntegrationRule * ApplyToKnotIntervals(KnotVector const &kv) const
Return an integration rule for KnotVector kv, defined by applying this rule on each knot interval...
Definition: intrules.cpp:181
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.
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:419
Direct sparse solver using UMFPACK.
Definition: solvers.hpp:1070
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition: gridfunc.cpp:366
void SetPatchRules1D(const int patch, std::vector< const IntegrationRule *> &ir1D)
Set 1D integration rules to be used as a tensor product rule on the patch with index patch...
Definition: intrules.cpp:2031
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:50
char vishost[]
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:302
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
Timing object.
Definition: tic_toc.hpp:35
virtual const char * Name() const
Definition: fe_coll.hpp:80
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, double RTOLERANCE, double ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition: solvers.cpp:898
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
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition: tic_toc.cpp:408
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
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
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
void GetPatchKnotVectors(int p, Array< KnotVector *> &kv)
Definition: nurbs.cpp:2670
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:3194
Wrapper for AlgebraicMultigrid object.
Definition: algebraic.hpp:185
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...
void SetNURBSPatchIntRule(NURBSMeshRules *pr)
For patchwise integration, SetNURBSPatchIntRule must be called.
Definition: nonlininteg.hpp:62
void DegreeElevate(int rel_degree, int degree=16)
Definition: mesh.cpp:5204
int dim
Definition: ex24.cpp:53
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:835
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
Vector data type.
Definition: vector.hpp:58
virtual void Print(std::ostream &os=mfem::out) const
Definition: mesh.hpp:2011
Vector with associated FE space and LinearFormIntegrators.
Definition: linearform.hpp:24
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:8302
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:3696
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition: device.hpp:121
void AssembleAndSolve(LinearForm &b, BilinearFormIntegrator *bfi, Array< int > const &ess_tdof_list, const bool pa, const bool algebraic_ceed, GridFunction &x)
int GetNP() const
Definition: nurbs.hpp:406
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:403
Class for defining different integration rules on each NURBS patch.
Definition: intrules.hpp:275
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition: solvers.cpp:3099