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