MFEM  v3.3.2
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
geom.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #include "fem.hpp"
13 
14 namespace mfem
15 {
16 
17 const char *Geometry::Name[NumGeom] =
18 { "Point", "Segment", "Triangle", "Square", "Tetrahedron", "Cube" };
19 
20 const double Geometry::Volume[NumGeom] =
21 { 1.0, 1.0, 0.5, 1.0, 1./6, 1.0 };
22 
24 {
25  // Vertices for Geometry::POINT
26  GeomVert[0] = NULL; // No vertices, dimension is 0
27 
28  // Vertices for Geometry::SEGMENT
29  GeomVert[1] = new IntegrationRule(2);
30 
31  GeomVert[1]->IntPoint(0).x = 0.0;
32  GeomVert[1]->IntPoint(1).x = 1.0;
33 
34  // Vertices for Geometry::TRIANGLE
35  GeomVert[2] = new IntegrationRule(3);
36 
37  GeomVert[2]->IntPoint(0).x = 0.0;
38  GeomVert[2]->IntPoint(0).y = 0.0;
39 
40  GeomVert[2]->IntPoint(1).x = 1.0;
41  GeomVert[2]->IntPoint(1).y = 0.0;
42 
43  GeomVert[2]->IntPoint(2).x = 0.0;
44  GeomVert[2]->IntPoint(2).y = 1.0;
45 
46  // Vertices for Geometry::SQUARE
47  GeomVert[3] = new IntegrationRule(4);
48 
49  GeomVert[3]->IntPoint(0).x = 0.0;
50  GeomVert[3]->IntPoint(0).y = 0.0;
51 
52  GeomVert[3]->IntPoint(1).x = 1.0;
53  GeomVert[3]->IntPoint(1).y = 0.0;
54 
55  GeomVert[3]->IntPoint(2).x = 1.0;
56  GeomVert[3]->IntPoint(2).y = 1.0;
57 
58  GeomVert[3]->IntPoint(3).x = 0.0;
59  GeomVert[3]->IntPoint(3).y = 1.0;
60 
61  // Vertices for Geometry::TETRAHEDRON
62  GeomVert[4] = new IntegrationRule(4);
63  GeomVert[4]->IntPoint(0).x = 0.0;
64  GeomVert[4]->IntPoint(0).y = 0.0;
65  GeomVert[4]->IntPoint(0).z = 0.0;
66 
67  GeomVert[4]->IntPoint(1).x = 1.0;
68  GeomVert[4]->IntPoint(1).y = 0.0;
69  GeomVert[4]->IntPoint(1).z = 0.0;
70 
71  GeomVert[4]->IntPoint(2).x = 0.0;
72  GeomVert[4]->IntPoint(2).y = 1.0;
73  GeomVert[4]->IntPoint(2).z = 0.0;
74 
75  GeomVert[4]->IntPoint(3).x = 0.0;
76  GeomVert[4]->IntPoint(3).y = 0.0;
77  GeomVert[4]->IntPoint(3).z = 1.0;
78 
79  // Vertices for Geometry::CUBE
80  GeomVert[5] = new IntegrationRule(8);
81 
82  GeomVert[5]->IntPoint(0).x = 0.0;
83  GeomVert[5]->IntPoint(0).y = 0.0;
84  GeomVert[5]->IntPoint(0).z = 0.0;
85 
86  GeomVert[5]->IntPoint(1).x = 1.0;
87  GeomVert[5]->IntPoint(1).y = 0.0;
88  GeomVert[5]->IntPoint(1).z = 0.0;
89 
90  GeomVert[5]->IntPoint(2).x = 1.0;
91  GeomVert[5]->IntPoint(2).y = 1.0;
92  GeomVert[5]->IntPoint(2).z = 0.0;
93 
94  GeomVert[5]->IntPoint(3).x = 0.0;
95  GeomVert[5]->IntPoint(3).y = 1.0;
96  GeomVert[5]->IntPoint(3).z = 0.0;
97 
98  GeomVert[5]->IntPoint(4).x = 0.0;
99  GeomVert[5]->IntPoint(4).y = 0.0;
100  GeomVert[5]->IntPoint(4).z = 1.0;
101 
102  GeomVert[5]->IntPoint(5).x = 1.0;
103  GeomVert[5]->IntPoint(5).y = 0.0;
104  GeomVert[5]->IntPoint(5).z = 1.0;
105 
106  GeomVert[5]->IntPoint(6).x = 1.0;
107  GeomVert[5]->IntPoint(6).y = 1.0;
108  GeomVert[5]->IntPoint(6).z = 1.0;
109 
110  GeomVert[5]->IntPoint(7).x = 0.0;
111  GeomVert[5]->IntPoint(7).y = 1.0;
112  GeomVert[5]->IntPoint(7).z = 1.0;
113 
114  GeomCenter[POINT].x = 0.0;
115  GeomCenter[POINT].y = 0.0;
116  GeomCenter[POINT].z = 0.0;
117 
118  GeomCenter[SEGMENT].x = 0.5;
119  GeomCenter[SEGMENT].y = 0.0;
120  GeomCenter[SEGMENT].z = 0.0;
121 
122  GeomCenter[TRIANGLE].x = 1.0 / 3.0;
123  GeomCenter[TRIANGLE].y = 1.0 / 3.0;
124  GeomCenter[TRIANGLE].z = 0.0;
125 
126  GeomCenter[SQUARE].x = 0.5;
127  GeomCenter[SQUARE].y = 0.5;
128  GeomCenter[SQUARE].z = 0.0;
129 
130  GeomCenter[TETRAHEDRON].x = 0.25;
131  GeomCenter[TETRAHEDRON].y = 0.25;
132  GeomCenter[TETRAHEDRON].z = 0.25;
133 
134  GeomCenter[CUBE].x = 0.5;
135  GeomCenter[CUBE].y = 0.5;
136  GeomCenter[CUBE].z = 0.5;
137 
138  GeomToPerfGeomJac[POINT] = NULL;
139  GeomToPerfGeomJac[SEGMENT] = new DenseMatrix(1);
140  GeomToPerfGeomJac[TRIANGLE] = new DenseMatrix(2);
141  GeomToPerfGeomJac[SQUARE] = new DenseMatrix(2);
142  GeomToPerfGeomJac[TETRAHEDRON] = new DenseMatrix(3);
143  GeomToPerfGeomJac[CUBE] = new DenseMatrix(3);
144 
145  PerfGeomToGeomJac[POINT] = NULL;
146  PerfGeomToGeomJac[SEGMENT] = NULL;
147  PerfGeomToGeomJac[TRIANGLE] = new DenseMatrix(2);
148  PerfGeomToGeomJac[SQUARE] = NULL;
149  PerfGeomToGeomJac[TETRAHEDRON] = new DenseMatrix(3);
150  PerfGeomToGeomJac[CUBE] = NULL;
151 
152  GeomToPerfGeomJac[SEGMENT]->Diag(1.0, 1);
153  {
154  Linear2DFiniteElement TriFE;
156  tri_T.SetFE(&TriFE);
158  tri_T.FinalizeTransformation();
159  tri_T.SetIntPoint(&GeomCenter[TRIANGLE]);
160  *GeomToPerfGeomJac[TRIANGLE] = tri_T.Jacobian();
161  CalcInverse(tri_T.Jacobian(), *PerfGeomToGeomJac[TRIANGLE]);
162  }
163  GeomToPerfGeomJac[SQUARE]->Diag(1.0, 2);
164  {
165  Linear3DFiniteElement TetFE;
167  tet_T.SetFE(&TetFE);
169  tet_T.FinalizeTransformation();
170  tet_T.SetIntPoint(&GeomCenter[TETRAHEDRON]);
171  *GeomToPerfGeomJac[TETRAHEDRON] = tet_T.Jacobian();
172  CalcInverse(tet_T.Jacobian(), *PerfGeomToGeomJac[TETRAHEDRON]);
173  }
174  GeomToPerfGeomJac[CUBE]->Diag(1.0, 3);
175 }
176 
178 {
179  for (int i = 0; i < NumGeom; i++)
180  {
181  delete PerfGeomToGeomJac[i];
182  delete GeomToPerfGeomJac[i];
183  delete GeomVert[i];
184  }
185 }
186 
188 {
189  switch (GeomType)
190  {
191  case Geometry::POINT: return GeomVert[0];
192  case Geometry::SEGMENT: return GeomVert[1];
193  case Geometry::TRIANGLE: return GeomVert[2];
194  case Geometry::SQUARE: return GeomVert[3];
195  case Geometry::TETRAHEDRON: return GeomVert[4];
196  case Geometry::CUBE: return GeomVert[5];
197  default:
198  mfem_error ("Geometry::GetVertices(...)");
199  }
200  // make some compilers happy.
201  return GeomVert[0];
202 }
203 
204 // static method
206 {
207  switch (GeomType)
208  {
209  case Geometry::POINT:
210  ip.x = 0.0;
211  break;
212  case Geometry::SEGMENT:
213  ip.x = double(rand()) / RAND_MAX;
214  break;
215  case Geometry::TRIANGLE:
216  ip.x = double(rand()) / RAND_MAX;
217  ip.y = double(rand()) / RAND_MAX;
218  if (ip.x + ip.y > 1.0)
219  {
220  ip.x = 1.0 - ip.x;
221  ip.y = 1.0 - ip.y;
222  }
223  break;
224  case Geometry::SQUARE:
225  ip.x = double(rand()) / RAND_MAX;
226  ip.y = double(rand()) / RAND_MAX;
227  break;
229  ip.x = double(rand()) / RAND_MAX;
230  ip.y = double(rand()) / RAND_MAX;
231  ip.z = double(rand()) / RAND_MAX;
232  // map to the triangular prism obtained by extruding the reference
233  // triangle in z direction
234  if (ip.x + ip.y > 1.0)
235  {
236  ip.x = 1.0 - ip.x;
237  ip.y = 1.0 - ip.y;
238  }
239  // split the prism into 3 parts: 1 is the reference tet, and the
240  // other two tets (as given below) are mapped to the reference tet
241  if (ip.x + ip.z > 1.0)
242  {
243  // tet with vertices: (0,0,1),(1,0,1),(0,1,1),(1,0,0)
244  ip.x = ip.x + ip.z - 1.0;
245  // ip.y = ip.y;
246  ip.z = 1.0 - ip.z;
247  // mapped to: (0,0,0),(1,0,0),(0,1,0),(0,0,1)
248  }
249  else if (ip.x + ip.y + ip.z > 1.0)
250  {
251  // tet with vertices: (0,1,1),(0,1,0),(0,0,1),(1,0,0)
252  double x = ip.x;
253  ip.x = 1.0 - x - ip.z;
254  ip.y = 1.0 - x - ip.y;
255  ip.z = x;
256  // mapped to: (0,0,0),(1,0,0),(0,1,0),(0,0,1)
257  }
258  break;
259  case Geometry::CUBE:
260  ip.x = double(rand()) / RAND_MAX;
261  ip.y = double(rand()) / RAND_MAX;
262  ip.z = double(rand()) / RAND_MAX;
263  break;
264  default:
265  MFEM_ABORT("Unknown type of reference element!");
266  }
267 }
268 
269 
270 namespace internal
271 {
272 
273 // Fuzzy equality operator with absolute tolerance eps.
274 inline bool NearlyEqual(double x, double y, double eps)
275 {
276  return std::abs(x-y) <= eps;
277 }
278 
279 // Fuzzy greater than comparison operator with absolute tolerance eps.
280 // Returns true when x is greater than y by at least eps.
281 inline bool FuzzyGT(double x, double y, double eps)
282 {
283  return (x > y + eps);
284 }
285 
286 // Fuzzy less than comparison operator with absolute tolerance eps.
287 // Returns true when x is less than y by at least eps.
288 inline bool FuzzyLT(double x, double y, double eps)
289 {
290  return (x < y - eps);
291 }
292 
293 }
294 
295 // static method
296 bool Geometry::CheckPoint(int GeomType, const IntegrationPoint &ip)
297 {
298  switch (GeomType)
299  {
300  case Geometry::POINT:
301  if (ip.x != 0.0) { return false; }
302  break;
303  case Geometry::SEGMENT:
304  if (ip.x < 0.0 || ip.x > 1.0) { return false; }
305  break;
306  case Geometry::TRIANGLE:
307  if (ip.x < 0.0 || ip.y < 0.0 || ip.x+ip.y > 1.0) { return false; }
308  break;
309  case Geometry::SQUARE:
310  if (ip.x < 0.0 || ip.x > 1.0 || ip.y < 0.0 || ip.y > 1.0)
311  { return false; }
312  break;
314  if (ip.x < 0.0 || ip.y < 0.0 || ip.z < 0.0 ||
315  ip.x+ip.y+ip.z > 1.0) { return false; }
316  break;
317  case Geometry::CUBE:
318  if (ip.x < 0.0 || ip.x > 1.0 || ip.y < 0.0 || ip.y > 1.0 ||
319  ip.z < 0.0 || ip.z > 1.0) { return false; }
320  break;
321  default:
322  MFEM_ABORT("Unknown type of reference element!");
323  }
324  return true;
325 }
326 
327 bool Geometry::CheckPoint(int GeomType, const IntegrationPoint &ip, double eps)
328 {
329  switch (GeomType)
330  {
331  case Geometry::POINT:
332  if (! internal::NearlyEqual(ip.x, 0.0, eps))
333  {
334  return false;
335  }
336  break;
337  case Geometry::SEGMENT:
338  if ( internal::FuzzyLT(ip.x, 0.0, eps)
339  || internal::FuzzyGT(ip.x, 1.0, eps) )
340  {
341  return false;
342  }
343  break;
344  case Geometry::TRIANGLE:
345  if ( internal::FuzzyLT(ip.x, 0.0, eps)
346  || internal::FuzzyLT(ip.y, 0.0, eps)
347  || internal::FuzzyGT(ip.x+ip.y, 1.0, eps) )
348  {
349  return false;
350  }
351  break;
352  case Geometry::SQUARE:
353  if ( internal::FuzzyLT(ip.x, 0.0, eps)
354  || internal::FuzzyGT(ip.x, 1.0, eps)
355  || internal::FuzzyLT(ip.y, 0.0, eps)
356  || internal::FuzzyGT(ip.y, 1.0, eps) )
357  {
358  return false;
359  }
360  break;
362  if ( internal::FuzzyLT(ip.x, 0.0, eps)
363  || internal::FuzzyLT(ip.y, 0.0, eps)
364  || internal::FuzzyLT(ip.z, 0.0, eps)
365  || internal::FuzzyGT(ip.x+ip.y+ip.z, 1.0, eps) )
366  {
367  return false;
368  }
369  break;
370  case Geometry::CUBE:
371  if ( internal::FuzzyLT(ip.x, 0.0, eps)
372  || internal::FuzzyGT(ip.x, 1.0, eps)
373  || internal::FuzzyLT(ip.y, 0.0, eps)
374  || internal::FuzzyGT(ip.y, 1.0, eps)
375  || internal::FuzzyLT(ip.z, 0.0, eps)
376  || internal::FuzzyGT(ip.z, 1.0, eps) )
377  {
378  return false;
379  }
380  break;
381  default:
382  MFEM_ABORT("Unknown type of reference element!");
383  }
384  return true;
385 }
386 
387 
388 namespace internal
389 {
390 
391 template <int N, int dim>
392 inline bool IntersectSegment(double lbeg[N], double lend[N],
393  IntegrationPoint &end)
394 {
395  double t = 1.0;
396  bool out = false;
397  for (int i = 0; i < N; i++)
398  {
399  lbeg[i] = std::max(lbeg[i], 0.0); // remove round-off
400  if (lend[i] < 0.0)
401  {
402  out = true;
403  t = std::min(t, lbeg[i]/(lbeg[i]-lend[i]));
404  }
405  }
406  if (out)
407  {
408  if (dim >= 1) { end.x = t*lend[0] + (1.0-t)*lbeg[0]; }
409  if (dim >= 2) { end.y = t*lend[1] + (1.0-t)*lbeg[1]; }
410  if (dim >= 3) { end.z = t*lend[2] + (1.0-t)*lbeg[2]; }
411  return false;
412  }
413  return true;
414 }
415 
416 inline bool ProjectTriangle(double &x, double &y)
417 {
418  if (x < 0.0)
419  {
420  x = 0.0;
421  if (y < 0.0) { y = 0.0; }
422  else if (y > 1.0) { y = 1.0; }
423  return false;
424  }
425  if (y < 0.0)
426  {
427  if (x > 1.0) { x = 1.0; }
428  y = 0.0;
429  return false;
430  }
431  const double l3 = 1.0-x-y;
432  if (l3 < 0.0)
433  {
434  if (y - x > 1.0) { x = 0.0; y = 1.0; }
435  else if (y - x < -1.0) { x = 1.0; y = 0.0; }
436  else { x += l3/2; y += l3/2; }
437  return false;
438  }
439  return true;
440 }
441 
442 }
443 
444 // static method
445 bool Geometry::ProjectPoint(int GeomType, const IntegrationPoint &beg,
446  IntegrationPoint &end)
447 {
448  switch (GeomType)
449  {
450  case Geometry::POINT:
451  {
452  if (end.x != 0.0) { end.x = 0.0; return false; }
453  break;
454  }
455  case Geometry::SEGMENT:
456  {
457  if (end.x < 0.0) { end.x = 0.0; return false; }
458  if (end.x > 1.0) { end.x = 1.0; return false; }
459  break;
460  }
461  case Geometry::TRIANGLE:
462  {
463  double lend[3] = { end.x, end.y, 1-end.x-end.y };
464  double lbeg[3] = { beg.x, beg.y, 1-beg.x-beg.y };
465  return internal::IntersectSegment<3,2>(lbeg, lend, end);
466  }
467  case Geometry::SQUARE:
468  {
469  double lend[4] = { end.x, end.y, 1-end.x, 1.0-end.y };
470  double lbeg[4] = { beg.x, beg.y, 1-beg.x, 1.0-beg.y };
471  return internal::IntersectSegment<4,2>(lbeg, lend, end);
472  }
474  {
475  double lend[4] = { end.x, end.y, end.z, 1.0-end.x-end.y-end.z };
476  double lbeg[4] = { beg.x, beg.y, beg.z, 1.0-beg.x-beg.y-beg.z };
477  return internal::IntersectSegment<4,3>(lbeg, lend, end);
478  }
479  case Geometry::CUBE:
480  {
481  double lend[6] = { end.x, end.y, end.z,
482  1.0-end.x, 1.0-end.y, 1.0-end.z
483  };
484  double lbeg[6] = { beg.x, beg.y, beg.z,
485  1.0-beg.x, 1.0-beg.y, 1.0-beg.z
486  };
487  return internal::IntersectSegment<6,3>(lbeg, lend, end);
488  }
489  default:
490  MFEM_ABORT("Unknown type of reference element!");
491  }
492  return true;
493 }
494 
495 // static method
497 {
498  // If ip is outside the element, replace it with the point on the boundary
499  // that is closest to the original ip and return false; otherwise, return
500  // true without changing ip.
501 
502  switch (GeomType)
503  {
504  case SEGMENT:
505  {
506  if (ip.x < 0.0) { ip.x = 0.0; return false; }
507  else if (ip.x > 1.0) { ip.x = 1.0; return false; }
508  return true;
509  }
510 
511  case TRIANGLE:
512  {
513  return internal::ProjectTriangle(ip.x, ip.y);
514  }
515 
516  case SQUARE:
517  {
518  bool in_x, in_y;
519  if (ip.x < 0.0) { in_x = false; ip.x = 0.0; }
520  else if (ip.x > 1.0) { in_x = false; ip.x = 1.0; }
521  else { in_x = true; }
522  if (ip.y < 0.0) { in_y = false; ip.y = 0.0; }
523  else if (ip.y > 1.0) { in_y = false; ip.y = 1.0; }
524  else { in_y = true; }
525  return in_x && in_y;
526  }
527 
528  case TETRAHEDRON:
529  {
530  if (ip.z < 0.0)
531  {
532  ip.z = 0.0;
533  internal::ProjectTriangle(ip.x, ip.y);
534  return false;
535  }
536  if (ip.y < 0.0)
537  {
538  ip.y = 0.0;
539  internal::ProjectTriangle(ip.x, ip.z);
540  return false;
541  }
542  if (ip.x < 0.0)
543  {
544  ip.x = 0.0;
545  internal::ProjectTriangle(ip.y, ip.z);
546  return false;
547  }
548  const double l4 = 1.0-ip.x-ip.y-ip.z;
549  if (l4 < 0.0)
550  {
551  const double l4_3 = l4/3;
552  ip.x += l4_3;
553  ip.y += l4_3;
554  internal::ProjectTriangle(ip.x, ip.y);
555  ip.z = 1.0-ip.x-ip.y;
556  return false;
557  }
558  return true;
559  }
560 
561  case CUBE:
562  {
563  bool in_x, in_y, in_z;
564  if (ip.x < 0.0) { in_x = false; ip.x = 0.0; }
565  else if (ip.x > 1.0) { in_x = false; ip.x = 1.0; }
566  else { in_x = true; }
567  if (ip.y < 0.0) { in_y = false; ip.y = 0.0; }
568  else if (ip.y > 1.0) { in_y = false; ip.y = 1.0; }
569  else { in_y = true; }
570  if (ip.z < 0.0) { in_z = false; ip.z = 0.0; }
571  else if (ip.z > 1.0) { in_z = false; ip.z = 1.0; }
572  else { in_z = true; }
573  return in_x && in_y && in_z;
574  }
575 
576  default:
577  MFEM_ABORT("Reference element type is not supported!");
578  }
579  return true;
580 }
581 
582 void Geometry::GetPerfPointMat(int GeomType, DenseMatrix &pm)
583 {
584  switch (GeomType)
585  {
586  case Geometry::SEGMENT:
587  {
588  pm.SetSize (1, 2);
589  pm(0,0) = 0.0;
590  pm(0,1) = 1.0;
591  }
592  break;
593 
594  case Geometry::TRIANGLE:
595  {
596  pm.SetSize (2, 3);
597  pm(0,0) = 0.0; pm(1,0) = 0.0;
598  pm(0,1) = 1.0; pm(1,1) = 0.0;
599  pm(0,2) = 0.5; pm(1,2) = 0.86602540378443864676;
600  }
601  break;
602 
603  case Geometry::SQUARE:
604  {
605  pm.SetSize (2, 4);
606  pm(0,0) = 0.0; pm(1,0) = 0.0;
607  pm(0,1) = 1.0; pm(1,1) = 0.0;
608  pm(0,2) = 1.0; pm(1,2) = 1.0;
609  pm(0,3) = 0.0; pm(1,3) = 1.0;
610  }
611  break;
612 
614  {
615  pm.SetSize (3, 4);
616  pm(0,0) = 0.0; pm(1,0) = 0.0; pm(2,0) = 0.0;
617  pm(0,1) = 1.0; pm(1,1) = 0.0; pm(2,1) = 0.0;
618  pm(0,2) = 0.5; pm(1,2) = 0.86602540378443864676; pm(2,2) = 0.0;
619  pm(0,3) = 0.5; pm(1,3) = 0.28867513459481288225;
620  pm(2,3) = 0.81649658092772603273;
621  }
622  break;
623 
624  case Geometry::CUBE:
625  {
626  pm.SetSize (3, 8);
627  pm(0,0) = 0.0; pm(1,0) = 0.0; pm(2,0) = 0.0;
628  pm(0,1) = 1.0; pm(1,1) = 0.0; pm(2,1) = 0.0;
629  pm(0,2) = 1.0; pm(1,2) = 1.0; pm(2,2) = 0.0;
630  pm(0,3) = 0.0; pm(1,3) = 1.0; pm(2,3) = 0.0;
631  pm(0,4) = 0.0; pm(1,4) = 0.0; pm(2,4) = 1.0;
632  pm(0,5) = 1.0; pm(1,5) = 0.0; pm(2,5) = 1.0;
633  pm(0,6) = 1.0; pm(1,6) = 1.0; pm(2,6) = 1.0;
634  pm(0,7) = 0.0; pm(1,7) = 1.0; pm(2,7) = 1.0;
635  }
636  break;
637 
638  default:
639  mfem_error ("Geometry::GetPerfPointMat (...)");
640  }
641 }
642 
643 void Geometry::JacToPerfJac(int GeomType, const DenseMatrix &J,
644  DenseMatrix &PJ) const
645 {
646  if (PerfGeomToGeomJac[GeomType])
647  {
648  Mult(J, *PerfGeomToGeomJac[GeomType], PJ);
649  }
650  else
651  {
652  PJ = J;
653  }
654 }
655 
656 const int Geometry::NumBdrArray[NumGeom] = { 0, 2, 3, 4, 4, 6 };
657 const int Geometry::Dimension[NumGeom] = { 0, 1, 2, 2, 3, 3 };
658 const int Geometry::NumVerts[NumGeom] = { 1, 2, 3, 4, 4, 8 };
659 const int Geometry::NumEdges[NumGeom] = { 0, 1, 3, 4, 6, 12 };
660 const int Geometry::NumFaces[NumGeom] = { 0, 0, 1, 1, 4, 6 };
661 
662 const int Geometry::
664 const int Geometry::
666 
667 const int Geometry::
668 Constants<Geometry::SEGMENT>::Edges[1][2] = { {0, 1} };
669 const int Geometry::
670 Constants<Geometry::SEGMENT>::Orient[2][2] = { {0, 1}, {1, 0} };
671 const int Geometry::
673 
674 const int Geometry::
675 Constants<Geometry::TRIANGLE>::Edges[3][2] = {{0, 1}, {1, 2}, {2, 0}};
676 const int Geometry::
678 const int Geometry::
679 Constants<Geometry::TRIANGLE>::VertToVert::J[3][2] = {{1, 0}, {2, -3}, {2, 1}};
680 const int Geometry::
681 Constants<Geometry::TRIANGLE>::FaceVert[1][3] = {{0, 1, 2}};
682 const int Geometry::
684 {
685  {0, 1, 2}, {1, 0, 2}, {2, 0, 1},
686  {2, 1, 0}, {1, 2, 0}, {0, 2, 1}
687 };
688 const int Geometry::
689 Constants<Geometry::TRIANGLE>::InvOrient[6] = { 0, 1, 4, 3, 2, 5 };
690 
691 const int Geometry::
692 Constants<Geometry::SQUARE>::Edges[4][2] = {{0, 1}, {1, 2}, {2, 3}, {3, 0}};
693 const int Geometry::
695 const int Geometry::
697 {{1, 0}, {3, -4}, {2, 1}, {3, 2}};
698 const int Geometry::
699 Constants<Geometry::SQUARE>::FaceVert[1][4] = {{0, 1, 2, 3}};
700 const int Geometry::
702 {
703  {0, 1, 2, 3}, {0, 3, 2, 1}, {1, 2, 3, 0}, {1, 0, 3, 2},
704  {2, 3, 0, 1}, {2, 1, 0, 3}, {3, 0, 1, 2}, {3, 2, 1, 0}
705 };
706 const int Geometry::
707 Constants<Geometry::SQUARE>::InvOrient[8] = { 0, 1, 6, 3, 4, 5, 2, 7 };
708 
709 const int Geometry::
711 {{0, 1}, {0, 2}, {0, 3}, {1, 2}, {1, 3}, {2, 3}};
712 const int Geometry::
714 {
716  Geometry::TRIANGLE, Geometry::TRIANGLE
717 };
718 const int Geometry::
720 {{1, 2, 3}, {0, 3, 2}, {0, 1, 3}, {0, 2, 1}};
721 const int Geometry::
723 const int Geometry::
725 {{1, 0}, {2, 1}, {3, 2}, {2, 3}, {3, 4}, {3, 5}};
726 
727 const int Geometry::
729 {
730  {0, 1}, {1, 2}, {3, 2}, {0, 3}, {4, 5}, {5, 6},
731  {7, 6}, {4, 7}, {0, 4}, {1, 5}, {2, 6}, {3, 7}
732 };
733 const int Geometry::
735 {
737  Geometry::SQUARE, Geometry::SQUARE, Geometry::SQUARE
738 };
739 const int Geometry::
741 {
742  {3, 2, 1, 0}, {0, 1, 5, 4}, {1, 2, 6, 5},
743  {2, 3, 7, 6}, {3, 0, 4, 7}, {4, 5, 6, 7}
744 };
745 const int Geometry::
746 Constants<Geometry::CUBE>::VertToVert::I[8] = {0, 3, 5, 7, 8, 10, 11, 12};
747 const int Geometry::
749 {
750  {1, 0}, {3, 3}, {4, 8}, // 0,1:0 0,3:3 0,4:8
751  {2, 1}, {5, 9}, // 1,2:1 1,5:9
752  {3,-3}, {6,10}, // 2,3:-3 2,6:10
753  {7,11}, // 3,7:11
754  {5, 4}, {7, 7}, // 4,5:4 4,7:7
755  {6, 5}, // 5,6:5
756  {7,-7} // 6,7:-7
757 };
758 
760 
761 
763 {
765 }
766 
768 {
769  for (int i = 0; i < Geometry::NumGeom; i++)
770  {
771  for (int j = 0; j < RGeom[i].Size(); j++) { delete RGeom[i][j]; }
772  for (int j = 0; j < IntPts[i].Size(); j++) { delete IntPts[i][j]; }
773  }
774 }
775 
776 RefinedGeometry *GeometryRefiner::FindInRGeom(int Geom, int Times, int ETimes,
777  int Type)
778 {
779  Array<RefinedGeometry *> &RGA = RGeom[Geom];
780  for (int i = 0; i < RGA.Size(); i++)
781  {
782  RefinedGeometry &RG = *RGA[i];
783  if (RG.Times == Times && RG.ETimes == ETimes && RG.Type == Type)
784  {
785  return &RG;
786  }
787  }
788  return NULL;
789 }
790 
791 IntegrationRule *GeometryRefiner::FindInIntPts(int Geom, int NPts)
792 {
793  Array<IntegrationRule *> &IPA = IntPts[Geom];
794  for (int i = 0; i < IPA.Size(); i++)
795  {
796  IntegrationRule &ir = *IPA[i];
797  if (ir.GetNPoints() == NPts) { return &ir; }
798  }
799  return NULL;
800 }
801 
802 RefinedGeometry * GeometryRefiner::Refine(int Geom, int Times, int ETimes)
803 {
804  int i, j, k, l;
805 
806  Times = std::max(Times, 1);
807  ETimes = std::max(ETimes, 1);
808  const double *cp = poly1d.GetPoints(Times, type);
809 
810  RefinedGeometry *RG = FindInRGeom(Geom, Times, ETimes, type);
811  if (RG) { return RG; }
812 
813  switch (Geom)
814  {
815  case Geometry::SEGMENT:
816  {
817  RG = new RefinedGeometry(Times+1, 2*Times, 0);
818  RG->Times = Times;
819  RG->ETimes = 0;
820  RG->Type = type;
821  for (i = 0; i <= Times; i++)
822  {
823  IntegrationPoint &ip = RG->RefPts.IntPoint(i);
824  ip.x = cp[i];
825  }
826  Array<int> &G = RG->RefGeoms;
827  for (i = 0; i < Times; i++)
828  {
829  G[2*i+0] = i;
830  G[2*i+1] = i+1;
831  }
832 
833  RGeom[Geometry::SEGMENT].Append(RG);
834  return RG;
835  }
836 
837  case Geometry::TRIANGLE:
838  {
839  RG = new RefinedGeometry((Times+1)*(Times+2)/2, 3*Times*Times,
840  3*Times*(ETimes+1), 3*Times);
841  RG->Times = Times;
842  RG->ETimes = ETimes;
843  RG->Type = type;
844  for (k = j = 0; j <= Times; j++)
845  for (i = 0; i <= Times-j; i++, k++)
846  {
847  IntegrationPoint &ip = RG->RefPts.IntPoint(k);
848  ip.x = cp[i]/(cp[i] + cp[j] + cp[Times-i-j]);
849  ip.y = cp[j]/(cp[i] + cp[j] + cp[Times-i-j]);
850  }
851  Array<int> &G = RG->RefGeoms;
852  for (l = k = j = 0; j < Times; j++, k++)
853  for (i = 0; i < Times-j; i++, k++)
854  {
855  G[l++] = k;
856  G[l++] = k+1;
857  G[l++] = k+Times-j+1;
858  if (i+j+1 < Times)
859  {
860  G[l++] = k+1;
861  G[l++] = k+Times-j+2;
862  G[l++] = k+Times-j+1;
863  }
864  }
865  Array<int> &E = RG->RefEdges;
866  int lb = 0, li = 2*RG->NumBdrEdges;
867  // horizontal edges
868  for (k = 0; k < Times; k += Times/ETimes)
869  {
870  int &lt = (k == 0) ? lb : li;
871  j = k*(Times+1)-((k-1)*k)/2;
872  for (i = 0; i < Times-k; i++)
873  {
874  E[lt++] = j; j++;
875  E[lt++] = j;
876  }
877  }
878  // diagonal edges
879  for (k = Times; k > 0; k -= Times/ETimes)
880  {
881  int &lt = (k == Times) ? lb : li;
882  j = k;
883  for (i = 0; i < k; i++)
884  {
885  E[lt++] = j; j += Times-i;
886  E[lt++] = j;
887  }
888  }
889  // vertical edges
890  for (k = 0; k < Times; k += Times/ETimes)
891  {
892  int &lt = (k == 0) ? lb : li;
893  j = k;
894  for (i = 0; i < Times-k; i++)
895  {
896  E[lt++] = j; j += Times-i+1;
897  E[lt++] = j;
898  }
899  }
900 
901  RGeom[Geometry::TRIANGLE].Append(RG);
902  return RG;
903  }
904 
905  case Geometry::SQUARE:
906  {
907  RG = new RefinedGeometry((Times+1)*(Times+1), 4*Times*Times,
908  4*(ETimes+1)*Times, 4*Times);
909  RG->Times = Times;
910  RG->ETimes = ETimes;
911  RG->Type = type;
912  for (k = j = 0; j <= Times; j++)
913  for (i = 0; i <= Times; i++, k++)
914  {
915  IntegrationPoint &ip = RG->RefPts.IntPoint(k);
916  ip.x = cp[i];
917  ip.y = cp[j];
918  }
919  Array<int> &G = RG->RefGeoms;
920  for (l = k = j = 0; j < Times; j++, k++)
921  for (i = 0; i < Times; i++, k++)
922  {
923  G[l++] = k;
924  G[l++] = k+1;
925  G[l++] = k+Times+2;
926  G[l++] = k+Times+1;
927  }
928  Array<int> &E = RG->RefEdges;
929  int lb = 0, li = 2*RG->NumBdrEdges;
930  // horizontal edges
931  for (k = 0; k <= Times; k += Times/ETimes)
932  {
933  int &lt = (k == 0 || k == Times) ? lb : li;
934  for (i = 0, j = k*(Times+1); i < Times; i++)
935  {
936  E[lt++] = j; j++;
937  E[lt++] = j;
938  }
939  }
940  // vertical edges (in right-to-left order)
941  for (k = Times; k >= 0; k -= Times/ETimes)
942  {
943  int &lt = (k == Times || k == 0) ? lb : li;
944  for (i = 0, j = k; i < Times; i++)
945  {
946  E[lt++] = j; j += Times+1;
947  E[lt++] = j;
948  }
949  }
950 
951  RGeom[Geometry::SQUARE].Append(RG);
952  return RG;
953  }
954 
955  case Geometry::CUBE:
956  {
957  RG = new RefinedGeometry ((Times+1)*(Times+1)*(Times+1),
958  8*Times*Times*Times, 0);
959  RG->Times = Times;
960  RG->ETimes = ETimes;
961  RG->Type = type;
962  for (l = k = 0; k <= Times; k++)
963  for (j = 0; j <= Times; j++)
964  for (i = 0; i <= Times; i++, l++)
965  {
966  IntegrationPoint &ip = RG->RefPts.IntPoint(l);
967  ip.x = cp[i];
968  ip.y = cp[j];
969  ip.z = cp[k];
970  }
971  Array<int> &G = RG->RefGeoms;
972  for (l = k = 0; k < Times; k++)
973  for (j = 0; j < Times; j++)
974  for (i = 0; i < Times; i++)
975  {
976  G[l++] = i+0 + (j+0 + (k+0) * (Times+1)) * (Times+1);
977  G[l++] = i+1 + (j+0 + (k+0) * (Times+1)) * (Times+1);
978  G[l++] = i+1 + (j+1 + (k+0) * (Times+1)) * (Times+1);
979  G[l++] = i+0 + (j+1 + (k+0) * (Times+1)) * (Times+1);
980  G[l++] = i+0 + (j+0 + (k+1) * (Times+1)) * (Times+1);
981  G[l++] = i+1 + (j+0 + (k+1) * (Times+1)) * (Times+1);
982  G[l++] = i+1 + (j+1 + (k+1) * (Times+1)) * (Times+1);
983  G[l++] = i+0 + (j+1 + (k+1) * (Times+1)) * (Times+1);
984  }
985 
986  RGeom[Geometry::CUBE].Append(RG);
987  return RG;
988  }
989 
991  {
992  // subdivide the tetrahedron with vertices
993  // (0,0,0), (0,0,1), (1,1,1), (0,1,1)
994 
995  // vertices: 0 <= i <= j <= k <= Times
996  // (3-combination with repetitions)
997  // number of vertices: (n+3)*(n+2)*(n+1)/6, n = Times
998 
999  // elements: the vertices are: v1=(i,j,k), v2=v1+u1, v3=v2+u2, v4=v3+u3
1000  // where 0 <= i <= j <= k <= n-1 and
1001  // u1,u2,u3 is a permutation of (1,0,0),(0,1,0),(0,0,1)
1002  // such that all v2,v3,v4 have non-decreasing components
1003  // number of elements: n^3
1004 
1005  const int n = Times;
1006  RG = new RefinedGeometry((n+3)*(n+2)*(n+1)/6, 4*n*n*n, 0);
1007  RG->Times = Times;
1008  RG->ETimes = ETimes;
1009  RG->Type = type;
1010  // enumerate and define the vertices
1011  Array<int> vi((n+1)*(n+1)*(n+1));
1012  vi = -1;
1013  int m = 0;
1014  for (k = 0; k <= n; k++)
1015  for (j = 0; j <= k; j++)
1016  for (i = 0; i <= j; i++)
1017  {
1018  IntegrationPoint &ip = RG->RefPts.IntPoint(m);
1019  // map the coordinates to the reference tetrahedron
1020  // (0,0,0) -> (0,0,0)
1021  // (0,0,1) -> (1,0,0)
1022  // (1,1,1) -> (0,1,0)
1023  // (0,1,1) -> (0,0,1)
1024  double w = cp[k-j] + cp[i] + cp[j-i] + cp[Times-k];
1025  ip.x = cp[k-j]/w;
1026  ip.y = cp[i]/w;
1027  ip.z = cp[j-i]/w;
1028  l = i + (j + k * (n+1)) * (n+1);
1029  vi[l] = m;
1030  m++;
1031  }
1032  if (m != (n+3)*(n+2)*(n+1)/6)
1033  {
1034  mfem_error("GeometryRefiner::Refine() for TETRAHEDRON #1");
1035  }
1036  // elements
1037  Array<int> &G = RG->RefGeoms;
1038  m = 0;
1039  for (k = 0; k < n; k++)
1040  for (j = 0; j <= k; j++)
1041  for (i = 0; i <= j; i++)
1042  {
1043  // the ordering of the vertices is chosen to ensure:
1044  // 1) correct orientation
1045  // 2) the x,y,z edges are in the set of edges
1046  // {(0,1),(2,3), (0,2),(1,3)}
1047  // (goal is to ensure that subsequent refinement using
1048  // this procedure preserves the six tetrahedral shapes)
1049 
1050  // zyx: (i,j,k)-(i,j,k+1)-(i+1,j+1,k+1)-(i,j+1,k+1)
1051  G[m++] = vi[i+0 + (j+0 + (k+0) * (n+1)) * (n+1)];
1052  G[m++] = vi[i+0 + (j+0 + (k+1) * (n+1)) * (n+1)];
1053  G[m++] = vi[i+1 + (j+1 + (k+1) * (n+1)) * (n+1)];
1054  G[m++] = vi[i+0 + (j+1 + (k+1) * (n+1)) * (n+1)];
1055  if (j < k)
1056  {
1057  // yzx: (i,j,k)-(i+1,j+1,k+1)-(i,j+1,k)-(i,j+1,k+1)
1058  G[m++] = vi[i+0 + (j+0 + (k+0) * (n+1)) * (n+1)];
1059  G[m++] = vi[i+1 + (j+1 + (k+1) * (n+1)) * (n+1)];
1060  G[m++] = vi[i+0 + (j+1 + (k+0) * (n+1)) * (n+1)];
1061  G[m++] = vi[i+0 + (j+1 + (k+1) * (n+1)) * (n+1)];
1062  // yxz: (i,j,k)-(i,j+1,k)-(i+1,j+1,k+1)-(i+1,j+1,k)
1063  G[m++] = vi[i+0 + (j+0 + (k+0) * (n+1)) * (n+1)];
1064  G[m++] = vi[i+0 + (j+1 + (k+0) * (n+1)) * (n+1)];
1065  G[m++] = vi[i+1 + (j+1 + (k+1) * (n+1)) * (n+1)];
1066  G[m++] = vi[i+1 + (j+1 + (k+0) * (n+1)) * (n+1)];
1067  }
1068  if (i < j)
1069  {
1070  // xzy: (i,j,k)-(i+1,j,k)-(i+1,j+1,k+1)-(i+1,j,k+1)
1071  G[m++] = vi[i+0 + (j+0 + (k+0) * (n+1)) * (n+1)];
1072  G[m++] = vi[i+1 + (j+0 + (k+0) * (n+1)) * (n+1)];
1073  G[m++] = vi[i+1 + (j+1 + (k+1) * (n+1)) * (n+1)];
1074  G[m++] = vi[i+1 + (j+0 + (k+1) * (n+1)) * (n+1)];
1075  if (j < k)
1076  {
1077  // xyz: (i,j,k)-(i+1,j+1,k+1)-(i+1,j,k)-(i+1,j+1,k)
1078  G[m++] = vi[i+0 + (j+0 + (k+0) * (n+1)) * (n+1)];
1079  G[m++] = vi[i+1 + (j+1 + (k+1) * (n+1)) * (n+1)];
1080  G[m++] = vi[i+1 + (j+0 + (k+0) * (n+1)) * (n+1)];
1081  G[m++] = vi[i+1 + (j+1 + (k+0) * (n+1)) * (n+1)];
1082  }
1083  // zxy: (i,j,k)-(i+1,j+1,k+1)-(i,j,k+1)-(i+1,j,k+1)
1084  G[m++] = vi[i+0 + (j+0 + (k+0) * (n+1)) * (n+1)];
1085  G[m++] = vi[i+1 + (j+1 + (k+1) * (n+1)) * (n+1)];
1086  G[m++] = vi[i+0 + (j+0 + (k+1) * (n+1)) * (n+1)];
1087  G[m++] = vi[i+1 + (j+0 + (k+1) * (n+1)) * (n+1)];
1088  }
1089  }
1090  if (m != 4*n*n*n)
1091  {
1092  mfem_error("GeometryRefiner::Refine() for TETRAHEDRON #2");
1093  }
1094  for (i = 0; i < m; i++)
1095  if (G[i] < 0)
1096  {
1097  mfem_error("GeometryRefiner::Refine() for TETRAHEDRON #3");
1098  }
1099 
1100  RGeom[Geometry::TETRAHEDRON].Append(RG);
1101  return RG;
1102  }
1103 
1104  default:
1105 
1106  return NULL;
1107  }
1108 }
1109 
1111 {
1112  IntegrationRule *ir = NULL;
1113 
1114  switch (Geom)
1115  {
1116  case Geometry::SEGMENT:
1117  {
1118  if (Times < 2)
1119  {
1120  return NULL;
1121  }
1122  ir = FindInIntPts(Geom, Times-1);
1123  if (ir == NULL)
1124  {
1125  ir = new IntegrationRule(Times-1);
1126  for (int i = 1; i < Times; i++)
1127  {
1128  IntegrationPoint &ip = ir->IntPoint(i-1);
1129  ip.x = double(i) / Times;
1130  ip.y = ip.z = 0.0;
1131  }
1132  }
1133  }
1134  break;
1135 
1136  case Geometry::TRIANGLE:
1137  {
1138  if (Times < 3)
1139  {
1140  return NULL;
1141  }
1142  ir = FindInIntPts(Geom, ((Times-1)*(Times-2))/2);
1143  if (ir == NULL)
1144  {
1145  ir = new IntegrationRule(((Times-1)*(Times-2))/2);
1146  for (int k = 0, j = 1; j < Times-1; j++)
1147  for (int i = 1; i < Times-j; i++, k++)
1148  {
1149  IntegrationPoint &ip = ir->IntPoint(k);
1150  ip.x = double(i) / Times;
1151  ip.y = double(j) / Times;
1152  ip.z = 0.0;
1153  }
1154  }
1155  }
1156  break;
1157 
1158  case Geometry::SQUARE:
1159  {
1160  if (Times < 2)
1161  {
1162  return NULL;
1163  }
1164  ir = FindInIntPts(Geom, (Times-1)*(Times-1));
1165  if (ir == NULL)
1166  {
1167  ir = new IntegrationRule((Times-1)*(Times-1));
1168  for (int k = 0, j = 1; j < Times; j++)
1169  for (int i = 1; i < Times; i++, k++)
1170  {
1171  IntegrationPoint &ip = ir->IntPoint(k);
1172  ip.x = double(i) / Times;
1173  ip.y = double(j) / Times;
1174  ip.z = 0.0;
1175  }
1176  }
1177  }
1178  break;
1179 
1180  default:
1181  mfem_error("GeometryRefiner::RefineInterior(...)");
1182  }
1183 
1184  if (ir) { IntPts[Geom].Append(ir); }
1185  return ir;
1186 }
1187 
1189 
1190 }
int Size() const
Logical size of the array.
Definition: array.hpp:110
static const int InvOrient[NumOrient]
Definition: geom.hpp:107
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:83
DenseMatrix & GetPointMat()
Read and write access to the underlying point matrix describing the transformation.
Definition: eltrans.hpp:312
static const int NumGeom
Definition: geom.hpp:34
void JacToPerfJac(int GeomType, const DenseMatrix &J, DenseMatrix &PJ) const
Definition: geom.cpp:643
static void GetRandomPoint(int GeomType, IntegrationPoint &ip)
Get a random point in the reference element specified by GeomType.
Definition: geom.cpp:205
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:468
static const int FaceVert[NumFaces][MaxFaceVert]
Definition: geom.hpp:177
void SetIntPoint(const IntegrationPoint *ip)
Definition: eltrans.hpp:52
RefinedGeometry * Refine(int Geom, int Times, int ETimes=1)
Definition: geom.cpp:802
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
static const double Volume[NumGeom]
Definition: geom.hpp:37
static const int NumEdges[NumGeom]
Definition: geom.hpp:40
Array< int > RefEdges
Definition: geom.hpp:211
const IntegrationRule * GetVertices(int GeomType)
Return an IntegrationRule consisting of all vertices of the given Geometry::Type, GeomType...
Definition: geom.cpp:187
const double * GetPoints(const int p, const int type)
Get the coordinates of the points of the given Quadrature1D type.
Definition: fe.cpp:6552
static const int NumFaces[NumGeom]
Definition: geom.hpp:41
static const int InvOrient[NumOrient]
Definition: geom.hpp:145
Geometry Geometries
Definition: geom.cpp:759
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:225
static const int Edges[NumEdges][2]
Definition: geom.hpp:127
int dim
Definition: ex3.cpp:47
static const int Dimension[NumGeom]
Definition: geom.hpp:38
const IntegrationRule * RefineInterior(int Geom, int Times)
Definition: geom.cpp:1110
static const int FaceTypes[NumFaces]
Definition: geom.hpp:193
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:67
static const int NumVerts[NumGeom]
Definition: geom.hpp:39
static const int Edges[NumEdges][2]
Definition: geom.hpp:115
GeometryRefiner GlobGeometryRefiner
Definition: geom.cpp:1188
Class for linear FE on tetrahedron.
Definition: fe.hpp:800
IntegrationRule RefPts
Definition: geom.hpp:210
static const int NumBdrArray[NumGeom]
Definition: geom.hpp:35
Class for linear FE on triangle.
Definition: fe.hpp:520
static const int InvOrient[NumOrient]
Definition: geom.hpp:165
static const int Orient[NumOrient][NumVert]
Definition: geom.hpp:118
static const int Edges[NumEdges][2]
Definition: geom.hpp:173
static const char * Name[NumGeom]
Definition: geom.hpp:36
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:3156
static const int InvOrient[NumOrient]
Definition: geom.hpp:119
static const int Edges[NumEdges][2]
Definition: geom.hpp:153
static bool ProjectPoint(int GeomType, const IntegrationPoint &beg, IntegrationPoint &end)
Project a point end, onto the given Geometry::Type, GeomType.
Definition: geom.cpp:445
void GetPerfPointMat(int GeomType, DenseMatrix &pm)
Definition: geom.cpp:582
void mfem_error(const char *msg)
Definition: error.cpp:107
static const int FaceVert[NumFaces][MaxFaceVert]
Definition: geom.hpp:195
static const int FaceVert[NumFaces][NumVert]
Definition: geom.hpp:135
Class for integration point with weight.
Definition: intrules.hpp:25
aka closed Newton-Cotes
Definition: intrules.hpp:267
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
Definition: densemat.cpp:2371
static const int FaceVert[NumFaces][NumVert]
Definition: geom.hpp:161
static const int Orient[NumOrient][NumVert]
Definition: geom.hpp:164
void SetFE(const FiniteElement *FE)
Definition: eltrans.hpp:298
static const int Edges[NumEdges][2]
Definition: geom.hpp:191
static const int Orient[NumOrient][NumVert]
Definition: geom.hpp:106
static bool CheckPoint(int GeomType, const IntegrationPoint &ip)
Check if the given point is inside the given reference element.
Definition: geom.cpp:296
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:86
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:64
Array< int > RefGeoms
Definition: geom.hpp:211
static const int Orient[NumOrient][NumVert]
Definition: geom.hpp:143
Poly_1D poly1d
Definition: fe.cpp:6619