MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
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
26using namespace std;
27using namespace mfem;
28
30 Array<int> const& ess_tdof_list, const bool pa,
31 const bool algebraic_ceed, GridFunction & x);
32
33int 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#ifdef MFEM_USE_SINGLE
156 MFEM_ABORT("Reduced integration is not supported in single precision.");
157#endif
158
159 di->SetIntegrationMode(NonlinearFormIntegrator::Mode::PATCHWISE_REDUCED);
160 }
161 else if (patchAssembly)
162 {
163 di->SetIntegrationMode(NonlinearFormIntegrator::Mode::PATCHWISE);
164 }
165
166 NURBSMeshRules *patchRule = nullptr;
167 if (order < 0)
168 {
169 if (ir_order == -1) { ir_order = 2*fec->GetOrder(); }
170 cout << "Using ir_order " << ir_order << endl;
171
172 patchRule = new NURBSMeshRules(mesh.NURBSext->GetNP(), dim);
173 // Loop over patches and set a different rule for each patch.
174 for (int p=0; p<mesh.NURBSext->GetNP(); ++p)
175 {
177 mesh.NURBSext->GetPatchKnotVectors(p, kv);
178
179 std::vector<const IntegrationRule*> ir1D(dim);
180 const IntegrationRule *ir = &IntRules.Get(Geometry::SEGMENT, ir_order);
181
182 // Construct 1D integration rules by applying the rule ir to each
183 // knot span.
184 for (int i=0; i<dim; ++i)
185 {
186 ir1D[i] = ir->ApplyToKnotIntervals(*kv[i]);
187 }
188
189 patchRule->SetPatchRules1D(p, ir1D);
190 } // loop (p) over patches
191
192 patchRule->Finalize(mesh);
193 di->SetNURBSPatchIntRule(patchRule);
194 }
195
196 // 10. Assemble and solve the linear system
197 cout << "Assembling system patch-wise and solving" << endl;
198 AssembleAndSolve(b, di, ess_tdof_list, pa, algebraic_ceed, x);
199
200 delete patchRule;
201
202 // 11. Save the refined mesh and the solution. This output can be viewed
203 // later using GLVis: "glvis -m refined.mesh -g sol.gf".
204 ofstream mesh_ofs("refined.mesh");
205 mesh_ofs.precision(8);
206 mesh.Print(mesh_ofs);
207 ofstream sol_ofs("sol.gf");
208 sol_ofs.precision(8);
209 x.Save(sol_ofs);
210
211 // 12. Send the solution by socket to a GLVis server.
212 if (visualization)
213 {
214 char vishost[] = "localhost";
215 int visport = 19916;
216 socketstream sol_sock(vishost, visport);
217 sol_sock.precision(8);
218 sol_sock << "solution\n" << mesh << x << flush;
219 }
220
221 // 13. Optionally assemble element-wise and solve the linear system, to
222 // compare timings and compute relative error.
223 if (compareToElementWise)
224 {
225 Vector x_pw, x_ew;
226 x.GetTrueDofs(x_pw);
227
228 cout << "Assembling system element-wise and solving" << endl;
230 // Element-wise partial assembly is not supported on NURBS meshes, so we
231 // pass pa = false here.
232 AssembleAndSolve(b, d, ess_tdof_list, false, algebraic_ceed, x);
233
234 x.GetTrueDofs(x_ew);
235
236 const real_t solNorm = x_ew.Norml2();
237 x_ew -= x_pw;
238
239 cout << "Element-wise solution norm " << solNorm << endl;
240 cout << "Relative error of patch-wise solution "
241 << x_ew.Norml2() / solNorm << endl;
242 }
243
244 // 14. Free the used memory.
245 if (delete_fec)
246 {
247 delete fec;
248 }
249
250 return 0;
251}
252
253// This function deletes bfi when the BilinearForm goes out of scope.
255 Array<int> const& ess_tdof_list, const bool pa,
256 const bool algebraic_ceed, GridFunction & x)
257{
258 FiniteElementSpace *fespace = b.FESpace();
259 BilinearForm a(fespace);
260 if (pa) { a.SetAssemblyLevel(AssemblyLevel::PARTIAL); }
261
262 a.AddDomainIntegrator(bfi); // Takes ownership of bfi
263
264 StopWatch sw;
265 sw.Start();
266
267 // Assemble the bilinear form and the corresponding linear system, applying
268 // any necessary transformations such as: eliminating boundary conditions,
269 // applying conforming constraints for non-conforming AMR, etc.
270 a.Assemble();
271
272 sw.Stop();
273
274 const real_t timeAssemble = sw.RealTime();
275
276 sw.Clear();
277 sw.Start();
278
279 OperatorPtr A;
280 Vector B, X;
281 a.FormLinearSystem(ess_tdof_list, x, b, A, X, B);
282
283 sw.Stop();
284
285 const real_t timeFormLinearSystem = sw.RealTime();
286
287 cout << "Timing for Assemble: " << timeAssemble << " seconds" << endl;
288 cout << "Timing for FormLinearSystem: " << timeFormLinearSystem << " seconds"
289 << endl;
290 cout << "Timing for entire setup: " << timeAssemble + timeFormLinearSystem
291 << " seconds" << endl;
292
293 sw.Clear();
294 sw.Start();
295
296 // Solve the linear system A X = B.
297 if (!pa)
298 {
299#ifndef MFEM_USE_SUITESPARSE
300 // Use a simple symmetric Gauss-Seidel preconditioner with PCG.
301 GSSmoother M((SparseMatrix&)(*A));
302 PCG(*A, M, B, X, 1, 200, 1e-20, 0.0);
303#else
304 // If MFEM was compiled with SuiteSparse, use UMFPACK to solve the system.
305 UMFPackSolver umf_solver;
306 umf_solver.Control[UMFPACK_ORDERING] = UMFPACK_ORDERING_METIS;
307 umf_solver.SetOperator(*A);
308 umf_solver.Mult(B, X);
309#endif
310 }
311 else
312 {
313 if (UsesTensorBasis(*fespace))
314 {
315 if (algebraic_ceed)
316 {
317 ceed::AlgebraicSolver M(a, ess_tdof_list);
318 PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
319 }
320 else
321 {
322 OperatorJacobiSmoother M(a, ess_tdof_list);
323 PCG(*A, M, B, X, 1, 400, 1e-12, 0.0);
324 }
325 }
326 else
327 {
328 CG(*A, B, X, 1, 400, 1e-20, 0.0);
329 }
330 }
331
332 sw.Stop();
333 cout << "Timing for solve " << sw.RealTime() << endl;
334
335 // Recover the solution as a finite element grid function.
336 a.RecoverFEMSolution(X, b, x);
337}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
int Size() const
Return the logical size of the array.
Definition array.hpp:144
Abstract base class BilinearFormIntegrator.
A "square matrix" operator for the associated FE space and BLFIntegrators The sum of all the BLFInteg...
A coefficient that is constant across space and time.
The MFEM Device class abstracts hardware devices such as GPUs, as well as programming models such as ...
Definition device.hpp:123
void Print(std::ostream &out=mfem::out)
Print the configuration of the MFEM virtual device object.
Definition device.cpp:286
Class for domain integration .
Definition lininteg.hpp:109
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition fe_coll.hpp:27
int GetOrder() const
Return the order (polynomial degree) of the FE collection, corresponding to the order/degree returned...
Definition fe_coll.hpp:225
virtual const char * Name() const
Definition fe_coll.hpp:79
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition fespace.hpp:716
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition fespace.cpp:588
Data type for Gauss-Seidel smoother of sparse matrix.
Class for grid function - Vector with associated FE space.
Definition gridfunc.hpp:31
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
void GetTrueDofs(Vector &tv) const
Extract the true-dofs from the GridFunction.
Definition gridfunc.cpp:366
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
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
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Vector with associated FE space and LinearFormIntegrators.
Mesh data type.
Definition mesh.hpp:56
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition mesh.hpp:290
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition mesh.hpp:2288
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
void GetNodes(Vector &node_coord) const
Definition mesh.cpp:8973
void DegreeElevate(int rel_degree, int degree=16)
Definition mesh.cpp:5779
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition mesh.cpp:10970
int GetNP() const
Definition nurbs.hpp:482
void GetPatchKnotVectors(int p, Array< KnotVector * > &kv)
Definition nurbs.cpp:3136
Class for defining different integration rules on each NURBS patch.
Definition intrules.hpp:279
void Finalize(Mesh const &mesh)
Finalize() must be called before this class can be used for assembly. In particular,...
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....
void SetNURBSPatchIntRule(NURBSMeshRules *pr)
For patchwise integration, SetNURBSPatchIntRule must be called.
Pointer to an Operator of a specified type.
Definition handle.hpp:34
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:313
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
void PrintUsage(std::ostream &out) const
Print the usage message.
void PrintOptions(std::ostream &out) const
Print the options.
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 'var' to receive the value. Enable/disable tags are used to set the bool...
Definition optparser.hpp:82
bool Good() const
Return true if the command line options were parsed successfully.
Data type sparse matrix.
Definition sparsemat.hpp:51
Timing object.
Definition tic_toc.hpp:36
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition tic_toc.cpp:432
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition tic_toc.cpp:411
void Stop()
Stop the stopwatch.
Definition tic_toc.cpp:422
void Clear()
Clear the elapsed time on the stopwatch and restart it if it's running.
Definition tic_toc.cpp:406
Direct sparse solver using UMFPACK.
Definition solvers.hpp:1096
virtual void SetOperator(const Operator &op)
Factorize the given Operator op which must be a SparseMatrix.
Definition solvers.cpp:3102
real_t Control[UMFPACK_CONTROL]
Definition solvers.hpp:1106
virtual void Mult(const Vector &b, Vector &x) const
Direct solution of the linear system using UMFPACK.
Definition solvers.cpp:3197
Vector data type.
Definition vector.hpp:80
real_t Norml2() const
Returns the l2 norm of the vector.
Definition vector.cpp:832
Wrapper for AlgebraicMultigrid object.
int dim
Definition ex24.cpp:53
int main()
real_t b
Definition lissajous.cpp:42
real_t a
Definition lissajous.cpp:41
const int visport
void PCG(const Operator &A, Solver &B, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Preconditioned conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:913
void CG(const Operator &A, const Vector &b, Vector &x, int print_iter, int max_num_iter, real_t RTOLERANCE, real_t ATOLERANCE)
Conjugate gradient method. (tolerances are squared)
Definition solvers.cpp:898
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Definition fespace.hpp:1345
float real_t
Definition config.hpp:43
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486
const char vishost[]
real_t p(const Vector &x, real_t t)
void AssembleAndSolve(LinearForm &b, BilinearFormIntegrator *bfi, Array< int > const &ess_tdof_list, const bool pa, const bool algebraic_ceed, GridFunction &x)