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