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