MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
polar-nc.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11//
12// ----------------------------------------------
13// Polar NC: Generate polar non-conforming meshes
14// ----------------------------------------------
15//
16// This miniapp generates a circular sector mesh that consist of quadrilaterals
17// and triangles of similar sizes. The 3D version of the mesh is made of prisms
18// and tetrahedra. The mesh is non-conforming by design, and can optionally be
19// made curvilinear. The elements are ordered along a space-filling curve by
20// default, which makes the mesh ready for parallel non-conforming AMR in MFEM.
21//
22// The implementation also demonstrates how to initialize a non-conforming mesh
23// on the fly by marking hanging nodes with Mesh::AddVertexParents.
24//
25// Compile with: make polar-nc
26//
27// Sample runs: polar-nc --radius 1 --nsteps 10
28// polar-nc --aspect 2
29// polar-nc --dim 3 --order 4
30
31#include "mfem.hpp"
32#include <fstream>
33#include <iostream>
34
35using namespace mfem;
36using namespace std;
37
38
39struct Params2
40{
41 real_t r, dr;
42 real_t a, da;
43
44 Params2() = default;
45 Params2(real_t r0, real_t r1, real_t a0, real_t a1)
46 : r(r0), dr(r1 - r0), a(a0), da(a1 - a0) {}
47};
48
49Mesh* Make2D(int nsteps, real_t rstep, real_t phi, real_t aspect, int order,
50 bool sfc)
51{
52 Mesh *mesh = new Mesh(2, 0, 0);
53
54 int origin = mesh->AddVertex(0.0, 0.0);
55
56 // n is the number of steps in the polar direction
57 int n = 1;
58 while (phi * rstep/2 / n * aspect > rstep) { n++; }
59
60 real_t r = rstep;
61 int first = mesh->AddVertex(r, 0.0);
62
63 Array<Params2> params;
65
66 // create triangles around the origin
67 real_t prev_alpha = 0.0;
68 for (int i = 0; i < n; i++)
69 {
70 real_t alpha = phi * (i+1) / n;
71 mesh->AddVertex(r*cos(alpha), r*sin(alpha));
72 mesh->AddTriangle(origin, first+i, first+i+1);
73
74 params.Append(Params2(0, r, prev_alpha, alpha));
75 prev_alpha = alpha;
76 }
77
78 mesh->AddBdrSegment(origin, first, 1);
79 mesh->AddBdrSegment(first+n, origin, 2);
80
81 for (int k = 1; k < nsteps; k++)
82 {
83 // m is the number of polar steps of the previous row
84 int m = n;
85 int prev_first = first;
86
87 real_t prev_r = r;
88 r += rstep;
89
90 if (phi * (r + prev_r)/2 / n * aspect < rstep * sqrt(2))
91 {
92 if (k == 1) { blocks.Append(Pair<int, int>(mesh->GetNE(), n)); }
93
94 first = mesh->AddVertex(r, 0.0);
95 mesh->AddBdrSegment(prev_first, first, 1);
96
97 // create a row of quads, same number as in previous row
98 prev_alpha = 0.0;
99 for (int i = 0; i < n; i++)
100 {
101 real_t alpha = phi * (i+1) / n;
102 mesh->AddVertex(r*cos(alpha), r*sin(alpha));
103 mesh->AddQuad(prev_first+i, first+i, first+i+1, prev_first+i+1);
104
105 params.Append(Params2(prev_r, r, prev_alpha, alpha));
106 prev_alpha = alpha;
107 }
108
109 mesh->AddBdrSegment(first+n, prev_first+n, 2);
110 }
111 else // we need to double the number of elements per row
112 {
113 n *= 2;
114
115 blocks.Append(Pair<int, int>(mesh->GetNE(), n));
116
117 // first create hanging vertices
118 int hang = 0; // init to suppress gcc warning
119 for (int i = 0; i < m; i++)
120 {
121 real_t alpha = phi * (2*i+1) / n;
122 int index = mesh->AddVertex(prev_r*cos(alpha), prev_r*sin(alpha));
123 mesh->AddVertexParents(index, prev_first+i, prev_first+i+1);
124 if (!i) { hang = index; }
125 }
126
127 first = mesh->AddVertex(r, 0.0);
128 int a = prev_first, b = first;
129
130 mesh->AddBdrSegment(a, b, 1);
131
132 // create a row of quad pairs
133 prev_alpha = 0.0;
134 for (int i = 0; i < m; i++)
135 {
136 int c = hang+i, e = a+1;
137
138 real_t alpha_half = phi * (2*i+1) / n;
139 int d = mesh->AddVertex(r*cos(alpha_half), r*sin(alpha_half));
140
141 real_t alpha = phi * (2*i+2) / n;
142 int f = mesh->AddVertex(r*cos(alpha), r*sin(alpha));
143
144 mesh->AddQuad(a, b, d, c);
145 mesh->AddQuad(c, d, f, e);
146
147 a = e, b = f;
148
149 params.Append(Params2(prev_r, r, prev_alpha, alpha_half));
150 params.Append(Params2(prev_r, r, alpha_half, alpha));
151 prev_alpha = alpha;
152 }
153
154 mesh->AddBdrSegment(b, a, 2);
155 }
156 }
157
158 for (int i = 0; i < n; i++)
159 {
160 mesh->AddBdrSegment(first+i, first+i+1, 3);
161 }
162
163 // reorder blocks of elements with Grid SFC ordering
164 if (sfc)
165 {
166 blocks.Append(Pair<int, int>(mesh->GetNE(), 0));
167
168 Array<Params2> new_params(params.Size());
169
170 Array<int> ordering(mesh->GetNE());
171 for (int i = 0; i < blocks[0].one; i++)
172 {
173 ordering[i] = i;
174 new_params[i] = params[i];
175 }
176
177 Array<int> coords;
178 for (int i = 0; i < blocks.Size()-1; i++)
179 {
180 int beg = blocks[i].one;
181 int width = blocks[i].two;
182 int height = (blocks[i+1].one - blocks[i].one) / width;
183
184 NCMesh::GridSfcOrdering2D(width, height, coords);
185
186 for (int j = 0, k = 0; j < coords.Size(); k++, j += 2)
187 {
188 int sfc_index = ((i & 1) ? coords[j] : (width-1 - coords[j]))
189 + coords[j+1]*width;
190 int old_index = beg + sfc_index;
191
192 ordering[old_index] = beg + k;
193 new_params[beg + k] = params[old_index];
194 }
195 }
196
197 mesh->ReorderElements(ordering, false);
198
199 mfem::Swap(params, new_params);
200 }
201
202 mesh->FinalizeMesh();
203
204 // create high-order curvature
205 if (order > 1)
206 {
207 mesh->SetCurvature(order);
208
209 GridFunction *nodes = mesh->GetNodes();
210 const FiniteElementSpace *fes = mesh->GetNodalFESpace();
211
212 Array<int> dofs;
213 MFEM_ASSERT(params.Size() == mesh->GetNE(), "");
214
215 for (int i = 0; i < mesh->GetNE(); i++)
216 {
217 const Params2 &par = params[i];
218 const IntegrationRule &ir = fes->GetFE(i)->GetNodes();
220 fes->GetElementDofs(i, dofs);
221
222 for (int j = 0; j < dofs.Size(); j++)
223 {
224 real_t a;
225 if (geom == Geometry::SQUARE)
226 {
227 r = par.r + ir[j].x * par.dr;
228 a = par.a + ir[j].y * par.da;
229 }
230 else
231 {
232 real_t rr = ir[j].x + ir[j].y;
233 if (std::abs(rr) < 1e-12) { continue; }
234 r = par.r + rr * par.dr;
235 a = par.a + ir[j].y/rr * par.da;
236 }
237 (*nodes)(fes->DofToVDof(dofs[j], 0)) = r*cos(a);
238 (*nodes)(fes->DofToVDof(dofs[j], 1)) = r*sin(a);
239 }
240 }
241
242 nodes->RestrictConforming();
243 }
244
245 return mesh;
246}
247
248
249const real_t pi2 = M_PI / 2;
250
251struct Params3
252{
253 real_t r, dr;
254 real_t u1, u2, u3;
255 real_t v1, v2, v3;
256
257 Params3() = default;
258 Params3(real_t r0, real_t r1,
259 real_t u1, real_t v1, real_t u2, real_t v2, real_t u3, real_t v3)
260 : r(r0), dr(r1 - r0), u1(u1), u2(u2), u3(u3), v1(v1), v2(v2), v3(v3) {}
261};
262
263struct Vert : public Hashed2
264{
265 int id;
266};
267
268int GetMidVertex(int v1, int v2, real_t r, real_t u, real_t v, bool hanging,
269 Mesh *mesh, HashTable<Vert> &hash)
270{
271 int vmid = hash.FindId(v1, v2);
272 if (vmid < 0)
273 {
274 vmid = hash.GetId(v1, v2);
275
276 real_t w = 1.0 - u - v;
277 real_t q = r / sqrt(u*u + v*v + w*w);
278 int index = mesh->AddVertex(u*q, v*q, w*q);
279
280 if (hanging) { mesh->AddVertexParents(index, v1, v2); }
281
282 hash[vmid].id = index;
283 }
284 return hash[vmid].id;
285}
286
287void MakeLayer(int vx1, int vy1, int vz1, int vx2, int vy2, int vz2, int level,
288 real_t r1, real_t r2, real_t u1, real_t v1, real_t u2, real_t v2,
289 real_t u3, real_t v3, bool bnd1, bool bnd2, bool bnd3, bool bnd4,
290 Mesh *mesh, HashTable<Vert> &hash, Array<Params3> &params)
291{
292 if (!level)
293 {
294 mesh->AddWedge(vx1, vy1, vz1, vx2, vy2, vz2);
295
296 if (bnd1) { mesh->AddBdrQuad(vx1, vy1, vy2, vx2, 1); }
297 if (bnd2) { mesh->AddBdrQuad(vy1, vz1, vz2, vy2, 2); }
298 if (bnd3) { mesh->AddBdrQuad(vz1, vx1, vx2, vz2, 3); }
299 if (bnd4) { mesh->AddBdrTriangle(vx2, vy2, vz2, 4); }
300
301 params.Append(Params3(r1, r2, u1, v1, u2, v2, u3, v3));
302 }
303 else
304 {
305 real_t u12 = (u1+u2)/2, v12 = (v1+v2)/2;
306 real_t u23 = (u2+u3)/2, v23 = (v2+v3)/2;
307 real_t u31 = (u3+u1)/2, v31 = (v3+v1)/2;
308
309 bool hang = (level == 1);
310
311 int vxy1 = GetMidVertex(vx1, vy1, r1, u12, v12, hang, mesh, hash);
312 int vyz1 = GetMidVertex(vy1, vz1, r1, u23, v23, hang, mesh, hash);
313 int vxz1 = GetMidVertex(vx1, vz1, r1, u31, v31, hang, mesh, hash);
314 int vxy2 = GetMidVertex(vx2, vy2, r2, u12, v12, false, mesh, hash);
315 int vyz2 = GetMidVertex(vy2, vz2, r2, u23, v23, false, mesh, hash);
316 int vxz2 = GetMidVertex(vx2, vz2, r2, u31, v31, false, mesh, hash);
317
318 MakeLayer(vx1, vxy1, vxz1, vx2, vxy2, vxz2, level-1,
319 r1, r2, u1, v1, u12, v12, u31, v31,
320 bnd1, false, bnd3, bnd4, mesh, hash, params);
321 MakeLayer(vxy1, vy1, vyz1, vxy2, vy2, vyz2, level-1,
322 r1, r2, u12, v12, u2, v2, u23, v23,
323 bnd1, bnd2, false, bnd4, mesh, hash, params);
324 MakeLayer(vxz1, vyz1, vz1, vxz2, vyz2, vz2, level-1,
325 r1, r2, u31, v31, u23, v23, u3, v3,
326 false, bnd2, bnd3, bnd4, mesh, hash, params);
327 MakeLayer(vyz1, vxz1, vxy1, vyz2, vxz2, vxy2, level-1,
328 r1, r2, u23, v23, u31, v31, u12, v12,
329 false, false, false, bnd4, mesh, hash, params);
330 }
331}
332
333void MakeCenter(int origin, int vx, int vy, int vz, int level, real_t r,
334 real_t u1, real_t v1, real_t u2, real_t v2, real_t u3, real_t v3,
335 bool bnd1, bool bnd2, bool bnd3, bool bnd4,
336 Mesh *mesh, HashTable<Vert> &hash, Array<Params3> &params)
337{
338 if (!level)
339 {
340 mesh->AddTet(origin, vx, vy, vz);
341
342 if (bnd1) { mesh->AddBdrTriangle(0, vy, vx, 1); }
343 if (bnd2) { mesh->AddBdrTriangle(0, vz, vy, 2); }
344 if (bnd3) { mesh->AddBdrTriangle(0, vx, vz, 3); }
345 if (bnd4) { mesh->AddBdrTriangle(vx, vy, vz, 4); }
346
347 params.Append(Params3(0, r, u1, v1, u2, v2, u3, v3));
348 }
349 else
350 {
351 real_t u12 = (u1+u2)/2, v12 = (v1+v2)/2;
352 real_t u23 = (u2+u3)/2, v23 = (v2+v3)/2;
353 real_t u31 = (u3+u1)/2, v31 = (v3+v1)/2;
354
355 int vxy = GetMidVertex(vx, vy, r, u12, v12, false, mesh, hash);
356 int vyz = GetMidVertex(vy, vz, r, u23, v23, false, mesh, hash);
357 int vxz = GetMidVertex(vx, vz, r, u31, v31, false, mesh, hash);
358
359 MakeCenter(origin, vx, vxy, vxz, level-1, r, u1, v1, u12, v12, u31, v31,
360 bnd1, false, bnd3, bnd4, mesh, hash, params);
361 MakeCenter(origin, vxy, vy, vyz, level-1, r, u12, v12, u2, v2, u23, v23,
362 bnd1, bnd2, false, bnd4, mesh, hash, params);
363 MakeCenter(origin, vxz, vyz, vz, level-1, r, u31, v31, u23, v23, u3, v3,
364 false, bnd2, bnd3, bnd4, mesh, hash, params);
365 MakeCenter(origin, vyz, vxz, vxy, level-1, r, u23, v23, u31, v31, u12, v12,
366 false, false, false, bnd4, mesh, hash, params);
367 }
368}
369
370Mesh* Make3D(int nsteps, real_t rstep, real_t aspect, int order, bool sfc)
371{
372 Mesh *mesh = new Mesh(3, 0, 0);
373
374 HashTable<Vert> hash;
375 Array<Params3> params;
376
377 int origin = mesh->AddVertex(0, 0, 0);
378
379 real_t r = rstep;
380 int a = mesh->AddVertex(r, 0, 0);
381 int b = mesh->AddVertex(0, r, 0);
382 int c = mesh->AddVertex(0, 0, r);
383
384 int levels = 0;
385 while (pi2 * rstep / (1 << levels) * aspect > rstep) { levels++; }
386
387 MakeCenter(origin, a, b, c, levels, r, 1, 0, 0, 1, 0, 0,
388 true, true, true, (nsteps == 1), mesh, hash, params);
389
390 for (int k = 1; k < nsteps; k++)
391 {
392 real_t prev_r = r;
393 r += rstep;
394
395 if ((prev_r + rstep/2) * pi2 * aspect / (1 << levels) > rstep * sqrt(2))
396 {
397 levels++;
398 }
399
400 int d = mesh->AddVertex(r, 0, 0);
401 int e = mesh->AddVertex(0, r, 0);
402 int f = mesh->AddVertex(0, 0, r);
403
404 MakeLayer(a, b, c, d, e, f, levels, prev_r, r,
405 1, 0, 0, 1, 0, 0, true, true, true, (k == nsteps-1),
406 mesh, hash, params);
407
408 a = d;
409 b = e;
410 c = f;
411 }
412
413 // reorder mesh with Hilbert spatial sort
414 if (sfc)
415 {
416 Array<int> ordering;
417 mesh->GetHilbertElementOrdering(ordering);
418 mesh->ReorderElements(ordering, false);
419
420 Array<Params3> new_params(params.Size());
421 for (int i = 0; i < ordering.Size(); i++)
422 {
423 new_params[ordering[i]] = params[i];
424 }
425 mfem::Swap(params, new_params);
426 }
427
428 mesh->FinalizeMesh();
429
430 // create high-order curvature
431 if (order > 1)
432 {
433 mesh->SetCurvature(order);
434
435 GridFunction *nodes = mesh->GetNodes();
436 const FiniteElementSpace *fes = mesh->GetNodalFESpace();
437
438 Array<int> dofs;
439 MFEM_ASSERT(params.Size() == mesh->GetNE(), "");
440
441 for (int i = 0; i < mesh->GetNE(); i++)
442 {
443 const Params3 &par = params[i];
444 const IntegrationRule &ir = fes->GetFE(i)->GetNodes();
446 fes->GetElementDofs(i, dofs);
447
448 for (int j = 0; j < dofs.Size(); j++)
449 {
450 const IntegrationPoint &ip = ir[j];
451
452 real_t u, v, w;
453 if (geom == Geometry::PRISM)
454 {
455 real_t l1 = 1.0 - ip.x - ip.y;
456 real_t l2 = ip.x, l3 = ip.y;
457 u = l1 * par.u1 + l2 * par.u2 + l3 * par.u3;
458 v = l1 * par.v1 + l2 * par.v2 + l3 * par.v3;
459 w = 1.0 - u - v;
460 r = par.r + ip.z * par.dr;
461 }
462 else
463 {
464 u = ip.x * par.u1 + ip.y * par.u2 + ip.z * par.u3;
465 v = ip.x * par.v1 + ip.y * par.v2 + ip.z * par.v3;
466 real_t rr = ip.x + ip.y + ip.z;
467 if (std::abs(rr) < 1e-12) { continue; }
468 w = rr - u - v;
469 r = par.r + rr * par.dr;
470 }
471
472 real_t q = r / sqrt(u*u + v*v + w*w);
473 (*nodes)(fes->DofToVDof(dofs[j], 0)) = u*q;
474 (*nodes)(fes->DofToVDof(dofs[j], 1)) = v*q;
475 (*nodes)(fes->DofToVDof(dofs[j], 2)) = w*q;
476 }
477 }
478
479 nodes->RestrictConforming();
480 }
481
482 return mesh;
483}
484
485
486int main(int argc, char *argv[])
487{
488 int dim = 2;
489 real_t radius = 1.0;
490 int nsteps = 10;
491 real_t angle = 90;
492 real_t aspect = 1.0;
493 int order = 2;
494 bool sfc = true;
495 bool visualization = true;
496
497 // parse command line
498 OptionsParser args(argc, argv);
499 args.AddOption(&dim, "-d", "--dim", "Mesh dimension (2 or 3).");
500 args.AddOption(&radius, "-r", "--radius", "Radius of the domain.");
501 args.AddOption(&nsteps, "-n", "--nsteps",
502 "Number of elements along the radial direction");
503 args.AddOption(&aspect, "-a", "--aspect",
504 "Target aspect ratio of the elements.");
505 args.AddOption(&angle, "-phi", "--phi", "Angular range (2D only).");
506 args.AddOption(&order, "-o", "--order",
507 "Polynomial degree of mesh curvature.");
508 args.AddOption(&sfc, "-sfc", "--sfc", "-no-sfc", "--no-sfc",
509 "Try to order elements along a space-filling curve.");
510 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
511 "--no-visualization",
512 "Enable or disable GLVis visualization.");
513 args.Parse();
514 if (!args.Good())
515 {
516 args.PrintUsage(cout);
517 return EXIT_FAILURE;
518 }
519 args.PrintOptions(cout);
520
521 // validate options
522 MFEM_VERIFY(radius > 0, "");
523 MFEM_VERIFY(aspect > 0, "");
524 MFEM_VERIFY(dim >= 2 && dim <= 3, "");
525 MFEM_VERIFY(angle > 0 && angle < 360, "");
526 MFEM_VERIFY(nsteps > 0, "");
527
528 real_t phi = angle * M_PI / 180;
529
530 // generate
531 Mesh *mesh;
532 if (dim == 2)
533 {
534 mesh = Make2D(nsteps, radius/nsteps, phi, aspect, order, sfc);
535 }
536 else
537 {
538 mesh = Make3D(nsteps, radius/nsteps, aspect, order, sfc);
539 }
540
541 // save the final mesh
542 ofstream ofs("polar-nc.mesh");
543 ofs.precision(8);
544 mesh->Print(ofs);
545
546 // output the mesh to GLVis
547 if (visualization)
548 {
549 char vishost[] = "localhost";
550 int visport = 19916;
551 socketstream sol_sock(vishost, visport);
552 sol_sock.precision(8);
553 sol_sock << "mesh\n" << *mesh << flush;
554 }
555
556 delete mesh;
557
558 return EXIT_SUCCESS;
559}
int Size() const
Return the logical size of the array.
Definition: array.hpp:144
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition: array.hpp:769
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:220
DofTransformation * GetElementDofs(int elem, Array< int > &dofs) const
Returns indices of degrees of freedom of element 'elem'. The returned indices are offsets into an ldo...
Definition: fespace.cpp:2846
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition: fespace.cpp:3168
int DofToVDof(int dof, int vd, int ndofs=-1) const
Compute a single vdof corresponding to the index dof and the vector index vd.
Definition: fespace.cpp:249
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:395
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:31
void RestrictConforming()
Definition: gridfunc.cpp:1907
int GetId(int p1, int p2)
Get id of item whose parents are p1, p2... Create it if it doesn't exist.
Definition: hash.hpp:611
int FindId(int p1, int p2) const
Find id of item whose parents are p1, p2... Return -1 if it doesn't exist.
Definition: hash.hpp:702
Class for integration point with weight.
Definition: intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:100
Mesh data type.
Definition: mesh.hpp:56
void FinalizeMesh(int refine=0, bool fix_orientation=true)
Finalize the construction of any type of Mesh.
Definition: mesh.cpp:3128
int AddBdrQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:2063
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Adds a quadrilateral to the mesh given by 4 vertices v1 through v4.
Definition: mesh.cpp:1743
int AddWedge(int v1, int v2, int v3, int v4, int v5, int v6, int attr=1)
Adds a wedge to the mesh given by 6 vertices v1 through v6.
Definition: mesh.cpp:1778
const FiniteElementSpace * GetNodalFESpace() const
Definition: mesh.cpp:6206
int AddTriangle(int v1, int v2, int v3, int attr=1)
Adds a triangle to the mesh given by 3 vertices v1 through v3.
Definition: mesh.cpp:1729
virtual void Print(std::ostream &os=mfem::out, const std::string &comments="") const
Definition: mesh.hpp:2288
int AddVertex(real_t x, real_t y=0.0, real_t z=0.0)
Definition: mesh.cpp:1658
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:1226
int AddTet(int v1, int v2, int v3, int v4, int attr=1)
Adds a tetrahedron to the mesh given by 4 vertices v1 through v4.
Definition: mesh.cpp:1757
void ReorderElements(const Array< int > &ordering, bool reorder_vertices=true)
Definition: mesh.cpp:2458
int AddBdrSegment(int v1, int v2, int attr=1)
Definition: mesh.cpp:2035
void AddVertexParents(int i, int p1, int p2)
Mark vertex i as nonconforming, with parent vertices p1 and p2.
Definition: mesh.cpp:1682
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:8973
int AddBdrTriangle(int v1, int v2, int v3, int attr=1)
Definition: mesh.cpp:2049
void GetHilbertElementOrdering(Array< int > &ordering)
Definition: mesh.cpp:2406
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Set the curvature of the mesh nodes using the given polynomial degree.
Definition: mesh.cpp:6211
Geometry::Type GetElementBaseGeometry(int i) const
Definition: mesh.hpp:1385
static void GridSfcOrdering2D(int width, int height, Array< int > &coords)
Definition: ncmesh.cpp:5062
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
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.
Definition: optparser.hpp:159
A pair of objects.
Definition: sort_pairs.hpp:24
const real_t radius
Definition: distance.cpp:107
const real_t alpha
Definition: ex15.cpp:369
int dim
Definition: ex24.cpp:53
int visport
char vishost[]
int main()
int index(int i, int j, int nx, int ny)
Definition: life.cpp:235
real_t b
Definition: lissajous.cpp:42
real_t a
Definition: lissajous.cpp:41
real_t f(const Vector &p)
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:648
real_t u(const Vector &xvec)
Definition: lor_mms.hpp:22
float real_t
Definition: config.hpp:43
int GetMidVertex(int v1, int v2, real_t r, real_t u, real_t v, bool hanging, Mesh *mesh, HashTable< Vert > &hash)
Definition: polar-nc.cpp:268
void MakeLayer(int vx1, int vy1, int vz1, int vx2, int vy2, int vz2, int level, real_t r1, real_t r2, real_t u1, real_t v1, real_t u2, real_t v2, real_t u3, real_t v3, bool bnd1, bool bnd2, bool bnd3, bool bnd4, Mesh *mesh, HashTable< Vert > &hash, Array< Params3 > &params)
Definition: polar-nc.cpp:287
Mesh * Make3D(int nsteps, real_t rstep, real_t aspect, int order, bool sfc)
Definition: polar-nc.cpp:370
const real_t pi2
Definition: polar-nc.cpp:249
Mesh * Make2D(int nsteps, real_t rstep, real_t phi, real_t aspect, int order, bool sfc)
Definition: polar-nc.cpp:49
void MakeCenter(int origin, int vx, int vy, int vz, int level, real_t r, real_t u1, real_t v1, real_t u2, real_t v2, real_t u3, real_t v3, bool bnd1, bool bnd2, bool bnd3, bool bnd4, Mesh *mesh, HashTable< Vert > &hash, Array< Params3 > &params)
Definition: polar-nc.cpp:333