libMesh
fe_hermite_shape_3D.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // C++ includes
20 
21 // Local includes
22 #include "libmesh/fe.h"
23 #include "libmesh/elem.h"
24 #include "libmesh/fe_interface.h"
25 #include "libmesh/number_lookups.h"
26 
27 namespace
28 {
29 using namespace libMesh;
30 
31 // Compute the static coefficients for an element
32 void hermite_compute_coefs(const Elem * elem, std::vector<std::vector<Real>> & dxdxi
33 
34 #ifdef DEBUG
35  , std::vector<Real> & dydxi, std::vector<Real> & dzdeta, std::vector<Real> & dxdzeta,
36  std::vector<Real> & dzdxi, std::vector<Real> & dxdeta, std::vector<Real> & dydzeta
37 #endif
38  )
39 {
40  const FEFamily mapping_family = FEMap::map_fe_type(*elem);
41  const Order mapping_order (elem->default_order());
42  const ElemType mapping_elem_type (elem->type());
43 
44  const FEType map_fe_type(mapping_order, mapping_family);
45 
46  const int n_mapping_shape_functions =
47  FEInterface::n_shape_functions (3, map_fe_type,
48  mapping_elem_type);
49 
50  static const Point dofpt[2] = {Point(-1,-1,-1), Point(1,1,1)};
51 
52  FEInterface::shape_deriv_ptr shape_deriv_ptr =
53  FEInterface::shape_deriv_function(3, map_fe_type);
54 
55  for (int p = 0; p != 2; ++p)
56  {
57  dxdxi[0][p] = 0;
58  dxdxi[1][p] = 0;
59  dxdxi[2][p] = 0;
60 #ifdef DEBUG
61  dydxi[p] = 0;
62  dzdeta[p] = 0;
63  dxdzeta[p] = 0;
64  dzdxi[p] = 0;
65  dxdeta[p] = 0;
66  dydzeta[p] = 0;
67 #endif
68  for (int i = 0; i != n_mapping_shape_functions; ++i)
69  {
70  const Real ddxi = shape_deriv_ptr
71  (elem, mapping_order, i, 0, dofpt[p], false);
72  const Real ddeta = shape_deriv_ptr
73  (elem, mapping_order, i, 1, dofpt[p], false);
74  const Real ddzeta = shape_deriv_ptr
75  (elem, mapping_order, i, 2, dofpt[p], false);
76 
77  // dxdeta, dxdzeta, dydxi, dydzeta, dzdxi, dzdeta should all
78  // be 0!
79  const Point & point_i = elem->point(i);
80  dxdxi[0][p] += point_i(0) * ddxi;
81  dxdxi[1][p] += point_i(1) * ddeta;
82  dxdxi[2][p] += point_i(2) * ddzeta;
83 #ifdef DEBUG
84  dydxi[p] += point_i(1) * ddxi;
85  dzdeta[p] += point_i(2) * ddeta;
86  dxdzeta[p] += point_i(0) * ddzeta;
87  dzdxi[p] += point_i(2) * ddxi;
88  dxdeta[p] += point_i(0) * ddeta;
89  dydzeta[p] += point_i(1) * ddzeta;
90 #endif
91  }
92 
93  // No singular elements!
94  libmesh_assert(dxdxi[0][p]);
95  libmesh_assert(dxdxi[1][p]);
96  libmesh_assert(dxdxi[2][p]);
97  // No non-rectilinear or non-axis-aligned elements!
98 #ifdef DEBUG
99  libmesh_assert_less (std::abs(dydxi[p]), TOLERANCE);
100  libmesh_assert_less (std::abs(dzdeta[p]), TOLERANCE);
101  libmesh_assert_less (std::abs(dxdzeta[p]), TOLERANCE);
102  libmesh_assert_less (std::abs(dzdxi[p]), TOLERANCE);
103  libmesh_assert_less (std::abs(dxdeta[p]), TOLERANCE);
104  libmesh_assert_less (std::abs(dydzeta[p]), TOLERANCE);
105 #endif
106  }
107 }
108 
109 
110 
111 Real hermite_bases_3D (std::vector<unsigned int> & bases1D,
112  const std::vector<std::vector<Real>> & dxdxi,
113  const Order & o,
114  unsigned int i)
115 {
116  bases1D.clear();
117  bases1D.resize(3,0);
118  Real coef = 1.0;
119 
120  unsigned int e = o-2;
121 
122  // Nodes
123  if (i < 64)
124  {
125  switch (i / 8)
126  {
127  case 0:
128  break;
129  case 1:
130  bases1D[0] = 1;
131  break;
132  case 2:
133  bases1D[0] = 1;
134  bases1D[1] = 1;
135  break;
136  case 3:
137  bases1D[1] = 1;
138  break;
139  case 4:
140  bases1D[2] = 1;
141  break;
142  case 5:
143  bases1D[0] = 1;
144  bases1D[2] = 1;
145  break;
146  case 6:
147  bases1D[0] = 1;
148  bases1D[1] = 1;
149  bases1D[2] = 1;
150  break;
151  case 7:
152  bases1D[1] = 1;
153  bases1D[2] = 1;
154  break;
155  default:
156  libmesh_error_msg("Invalid basis node = " << i/8);
157  }
158 
159  unsigned int basisnum = i%8;
160  switch (basisnum) // DoF type
161  {
162  case 0: // DoF = value at node
163  coef = 1.0;
164  break;
165  case 1: // DoF = x derivative at node
166  coef = dxdxi[0][bases1D[0]];
167  bases1D[0] += 2; break;
168  case 2: // DoF = y derivative at node
169  coef = dxdxi[1][bases1D[1]];
170  bases1D[1] += 2; break;
171  case 3: // DoF = xy derivative at node
172  coef = dxdxi[0][bases1D[0]] * dxdxi[1][bases1D[1]];
173  bases1D[0] += 2; bases1D[1] += 2; break;
174  case 4: // DoF = z derivative at node
175  coef = dxdxi[2][bases1D[2]];
176  bases1D[2] += 2; break;
177  case 5: // DoF = xz derivative at node
178  coef = dxdxi[0][bases1D[0]] * dxdxi[2][bases1D[2]];
179  bases1D[0] += 2; bases1D[2] += 2; break;
180  case 6: // DoF = yz derivative at node
181  coef = dxdxi[1][bases1D[1]] * dxdxi[2][bases1D[2]];
182  bases1D[1] += 2; bases1D[2] += 2; break;
183  case 7: // DoF = xyz derivative at node
184  coef = dxdxi[0][bases1D[0]] * dxdxi[1][bases1D[1]] * dxdxi[2][bases1D[2]];
185  bases1D[0] += 2; bases1D[1] += 2; bases1D[2] += 2; break;
186  default:
187  libmesh_error_msg("Invalid basisnum = " << basisnum);
188  }
189  }
190  // Edges
191  else if (i < 64 + 12*4*e)
192  {
193  unsigned int basisnum = (i - 64) % (4*e);
194  switch ((i - 64) / (4*e))
195  {
196  case 0:
197  bases1D[0] = basisnum / 4 + 4;
198  bases1D[1] = basisnum % 4 / 2 * 2;
199  bases1D[2] = basisnum % 2 * 2;
200  if (basisnum % 4 / 2)
201  coef *= dxdxi[1][0];
202  if (basisnum % 2)
203  coef *= dxdxi[2][0];
204  break;
205  case 1:
206  bases1D[0] = basisnum % 4 / 2 * 2 + 1;
207  bases1D[1] = basisnum / 4 + 4;
208  bases1D[2] = basisnum % 2 * 2;
209  if (basisnum % 4 / 2)
210  coef *= dxdxi[0][1];
211  if (basisnum % 2)
212  coef *= dxdxi[2][0];
213  break;
214  case 2:
215  bases1D[0] = basisnum / 4 + 4;
216  bases1D[1] = basisnum % 4 / 2 * 2 + 1;
217  bases1D[2] = basisnum % 2 * 2;
218  if (basisnum % 4 / 2)
219  coef *= dxdxi[1][1];
220  if (basisnum % 2)
221  coef *= dxdxi[2][0];
222  break;
223  case 3:
224  bases1D[0] = basisnum % 4 / 2 * 2;
225  bases1D[1] = basisnum / 4 + 4;
226  bases1D[2] = basisnum % 2 * 2;
227  if (basisnum % 4 / 2)
228  coef *= dxdxi[0][0];
229  if (basisnum % 2)
230  coef *= dxdxi[2][0];
231  break;
232  case 4:
233  bases1D[0] = basisnum % 4 / 2 * 2;
234  bases1D[1] = basisnum % 2 * 2;
235  bases1D[2] = basisnum / 4 + 4;
236  if (basisnum % 4 / 2)
237  coef *= dxdxi[0][0];
238  if (basisnum % 2)
239  coef *= dxdxi[1][0];
240  break;
241  case 5:
242  bases1D[0] = basisnum % 4 / 2 * 2 + 1;
243  bases1D[1] = basisnum % 2 * 2;
244  bases1D[2] = basisnum / 4 + 4;
245  if (basisnum % 4 / 2)
246  coef *= dxdxi[0][1];
247  if (basisnum % 2)
248  coef *= dxdxi[1][0];
249  break;
250  case 6:
251  bases1D[0] = basisnum % 4 / 2 * 2 + 1;
252  bases1D[1] = basisnum % 2 * 2 + 1;
253  bases1D[2] = basisnum / 4 + 4;
254  if (basisnum % 4 / 2)
255  coef *= dxdxi[0][1];
256  if (basisnum % 2)
257  coef *= dxdxi[1][1];
258  break;
259  case 7:
260  bases1D[0] = basisnum % 4 / 2 * 2;
261  bases1D[1] = basisnum % 2 * 2 + 1;
262  bases1D[2] = basisnum / 4 + 4;
263  if (basisnum % 4 / 2)
264  coef *= dxdxi[0][0];
265  if (basisnum % 2)
266  coef *= dxdxi[1][1];
267  break;
268  case 8:
269  bases1D[0] = basisnum / 4 + 4;
270  bases1D[1] = basisnum % 4 / 2 * 2;
271  bases1D[2] = basisnum % 2 * 2 + 1;
272  if (basisnum % 4 / 2)
273  coef *= dxdxi[1][0];
274  if (basisnum % 2)
275  coef *= dxdxi[2][1];
276  break;
277  case 9:
278  bases1D[0] = basisnum % 4 / 2 * 2 + 1;
279  bases1D[1] = basisnum / 4 + 4;
280  bases1D[2] = basisnum % 2 * 2;
281  if (basisnum % 4 / 2)
282  coef *= dxdxi[0][1];
283  if (basisnum % 2)
284  coef *= dxdxi[2][1];
285  break;
286  case 10:
287  bases1D[0] = basisnum / 4 + 4;
288  bases1D[1] = basisnum % 4 / 2 * 2 + 1;
289  bases1D[2] = basisnum % 2 * 2 + 1;
290  if (basisnum % 4 / 2)
291  coef *= dxdxi[1][1];
292  if (basisnum % 2)
293  coef *= dxdxi[2][1];
294  break;
295  case 11:
296  bases1D[0] = basisnum % 4 / 2 * 2;
297  bases1D[1] = basisnum / 4 + 4;
298  bases1D[2] = basisnum % 2 * 2 + 1;
299  if (basisnum % 4 / 2)
300  coef *= dxdxi[0][0];
301  if (basisnum % 2)
302  coef *= dxdxi[2][1];
303  break;
304  default:
305  libmesh_error_msg("Invalid basis node = " << (i - 64) / (4*e));
306  }
307  }
308  // Faces
309  else if (i < 64 + 12*4*e + 6*2*e*e)
310  {
311  unsigned int basisnum = (i - 64 - 12*4*e) % (2*e*e);
312  switch ((i - 64 - 12*4*e) / (2*e*e))
313  {
314  case 0:
315  bases1D[0] = square_number_column[basisnum / 2];
316  bases1D[1] = square_number_row[basisnum / 2];
317  bases1D[2] = basisnum % 2 * 2;
318  if (basisnum % 2)
319  coef *= dxdxi[2][0];
320  break;
321  case 1:
322  bases1D[0] = square_number_column[basisnum / 2];
323  bases1D[1] = basisnum % 2 * 2;
324  bases1D[2] = square_number_row[basisnum / 2];
325  if (basisnum % 2)
326  coef *= dxdxi[1][0];
327  break;
328  case 2:
329  bases1D[0] = basisnum % 2 * 2 + 1;
330  bases1D[1] = square_number_column[basisnum / 2];
331  bases1D[2] = square_number_row[basisnum / 2];
332  if (basisnum % 2)
333  coef *= dxdxi[0][1];
334  break;
335  case 3:
336  bases1D[0] = square_number_column[basisnum / 2];
337  bases1D[1] = basisnum % 2 * 2 + 1;
338  bases1D[2] = square_number_row[basisnum / 2];
339  if (basisnum % 2)
340  coef *= dxdxi[1][1];
341  break;
342  case 4:
343  bases1D[0] = basisnum % 2 * 2;
344  bases1D[1] = square_number_column[basisnum / 2];
345  bases1D[2] = square_number_row[basisnum / 2];
346  if (basisnum % 2)
347  coef *= dxdxi[0][0];
348  break;
349  case 5:
350  bases1D[0] = square_number_column[basisnum / 2];
351  bases1D[1] = square_number_row[basisnum / 2];
352  bases1D[2] = basisnum % 2 * 2 + 1;
353  if (basisnum % 2)
354  coef *= dxdxi[2][1];
355  break;
356  default:
357  libmesh_error_msg("Invalid basis node = " << (i - 64 - 12*4*e) / (2*e*e));
358  }
359  }
360  // Interior
361  else
362  {
363  unsigned int basisnum = i - 64 - 12*4*e;
364  bases1D[0] = cube_number_column[basisnum] + 4;
365  bases1D[1] = cube_number_row[basisnum] + 4;
366  bases1D[2] = cube_number_page[basisnum] + 4;
367  }
368 
369  // No singular elements
370  libmesh_assert(coef);
371  return coef;
372 }
373 
374 
375 } // end anonymous namespace
376 
377 
378 namespace libMesh
379 {
380 
381 
382 template <>
384  const Order,
385  const unsigned int,
386  const Point &)
387 {
388  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
389  return 0.;
390 }
391 
392 
393 
394 template <>
396  const Order order,
397  const unsigned int i,
398  const Point & p,
399  const bool add_p_level)
400 {
401  libmesh_assert(elem);
402 
403  std::vector<std::vector<Real>> dxdxi(3, std::vector<Real>(2, 0));
404 
405 #ifdef DEBUG
406  std::vector<Real> dydxi(2), dzdeta(2), dxdzeta(2);
407  std::vector<Real> dzdxi(2), dxdeta(2), dydzeta(2);
408 #endif //DEBUG
409 
410  hermite_compute_coefs(elem, dxdxi
411 #ifdef DEBUG
412  , dydxi, dzdeta, dxdzeta, dzdxi, dxdeta, dydzeta
413 #endif
414  );
415 
416  const ElemType type = elem->type();
417 
418  const Order totalorder =
419  static_cast<Order>(order + add_p_level * elem->p_level());
420 
421  switch (totalorder)
422  {
423  // 3rd-order tricubic Hermite functions
424  case THIRD:
425  {
426  switch (type)
427  {
428  case HEX8:
429  case HEX20:
430  case HEX27:
431  {
432  libmesh_assert_less (i, 64);
433 
434  std::vector<unsigned int> bases1D;
435 
436  Real coef = hermite_bases_3D(bases1D, dxdxi, totalorder, i);
437 
438  return coef *
439  FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
440  FEHermite<1>::hermite_raw_shape(bases1D[1],p(1)) *
441  FEHermite<1>::hermite_raw_shape(bases1D[2],p(2));
442  }
443  default:
444  libmesh_error_msg("ERROR: Unsupported element type " << type);
445  }
446  }
447  // by default throw an error
448  default:
449  libmesh_error_msg("ERROR: Unsupported polynomial order " << totalorder);
450  }
451 }
452 
453 
454 
455 template <>
457  const Order,
458  const unsigned int,
459  const unsigned int,
460  const Point &)
461 {
462  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
463  return 0.;
464 }
465 
466 
467 
468 template <>
470  const Order order,
471  const unsigned int i,
472  const unsigned int j,
473  const Point & p,
474  const bool add_p_level)
475 {
476  libmesh_assert(elem);
477  libmesh_assert (j == 0 || j == 1 || j == 2);
478 
479  std::vector<std::vector<Real>> dxdxi(3, std::vector<Real>(2, 0));
480 
481 #ifdef DEBUG
482  std::vector<Real> dydxi(2), dzdeta(2), dxdzeta(2);
483  std::vector<Real> dzdxi(2), dxdeta(2), dydzeta(2);
484 #endif //DEBUG
485 
486  hermite_compute_coefs(elem, dxdxi
487 #ifdef DEBUG
488  , dydxi, dzdeta, dxdzeta, dzdxi, dxdeta, dydzeta
489 #endif
490  );
491 
492  const ElemType type = elem->type();
493 
494  const Order totalorder =
495  static_cast<Order>(order + add_p_level * elem->p_level());
496 
497  switch (totalorder)
498  {
499  // 3rd-order tricubic Hermite functions
500  case THIRD:
501  {
502  switch (type)
503  {
504  case HEX8:
505  case HEX20:
506  case HEX27:
507  {
508  libmesh_assert_less (i, 64);
509 
510  std::vector<unsigned int> bases1D;
511 
512  Real coef = hermite_bases_3D(bases1D, dxdxi, totalorder, i);
513 
514  switch (j) // Derivative type
515  {
516  case 0:
517  return coef *
518  FEHermite<1>::hermite_raw_shape_deriv(bases1D[0],p(0)) *
519  FEHermite<1>::hermite_raw_shape(bases1D[1],p(1)) *
520  FEHermite<1>::hermite_raw_shape(bases1D[2],p(2));
521  break;
522  case 1:
523  return coef *
524  FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
525  FEHermite<1>::hermite_raw_shape_deriv(bases1D[1],p(1)) *
526  FEHermite<1>::hermite_raw_shape(bases1D[2],p(2));
527  break;
528  case 2:
529  return coef *
530  FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
531  FEHermite<1>::hermite_raw_shape(bases1D[1],p(1)) *
532  FEHermite<1>::hermite_raw_shape_deriv(bases1D[2],p(2));
533  break;
534  default:
535  libmesh_error_msg("Invalid shape function derivative j = " << j);
536  }
537 
538  }
539  default:
540  libmesh_error_msg("ERROR: Unsupported element type " << type);
541  }
542  }
543  // by default throw an error
544  default:
545  libmesh_error_msg("ERROR: Unsupported polynomial order " << totalorder);
546  }
547 }
548 
549 
550 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
551 
552 template <>
554  const Order,
555  const unsigned int,
556  const unsigned int,
557  const Point &)
558 {
559  libmesh_error_msg("Hermite elements require the real element \nto construct gradient-based degrees of freedom.");
560  return 0.;
561 }
562 
563 
564 template <>
566  const Order order,
567  const unsigned int i,
568  const unsigned int j,
569  const Point & p,
570  const bool add_p_level)
571 {
572  libmesh_assert(elem);
573 
574  std::vector<std::vector<Real>> dxdxi(3, std::vector<Real>(2, 0));
575 
576 #ifdef DEBUG
577  std::vector<Real> dydxi(2), dzdeta(2), dxdzeta(2);
578  std::vector<Real> dzdxi(2), dxdeta(2), dydzeta(2);
579 #endif //DEBUG
580 
581  hermite_compute_coefs(elem, dxdxi
582 #ifdef DEBUG
583  , dydxi, dzdeta, dxdzeta, dzdxi, dxdeta, dydzeta
584 #endif
585  );
586 
587  const ElemType type = elem->type();
588 
589  const Order totalorder =
590  static_cast<Order>(order + add_p_level * elem->p_level());
591 
592  switch (totalorder)
593  {
594  // 3rd-order tricubic Hermite functions
595  case THIRD:
596  {
597  switch (type)
598  {
599  case HEX8:
600  case HEX20:
601  case HEX27:
602  {
603  libmesh_assert_less (i, 64);
604 
605  std::vector<unsigned int> bases1D;
606 
607  Real coef = hermite_bases_3D(bases1D, dxdxi, totalorder, i);
608 
609  switch (j) // Derivative type
610  {
611  case 0:
612  return coef *
614  FEHermite<1>::hermite_raw_shape(bases1D[1],p(1)) *
615  FEHermite<1>::hermite_raw_shape(bases1D[2],p(2));
616  break;
617  case 1:
618  return coef *
619  FEHermite<1>::hermite_raw_shape_deriv(bases1D[0],p(0)) *
620  FEHermite<1>::hermite_raw_shape_deriv(bases1D[1],p(1)) *
621  FEHermite<1>::hermite_raw_shape(bases1D[2],p(2));
622  break;
623  case 2:
624  return coef *
625  FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
627  FEHermite<1>::hermite_raw_shape(bases1D[2],p(2));
628  break;
629  case 3:
630  return coef *
631  FEHermite<1>::hermite_raw_shape_deriv(bases1D[0],p(0)) *
632  FEHermite<1>::hermite_raw_shape(bases1D[1],p(1)) *
633  FEHermite<1>::hermite_raw_shape_deriv(bases1D[2],p(2));
634  break;
635  case 4:
636  return coef *
637  FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
638  FEHermite<1>::hermite_raw_shape_deriv(bases1D[1],p(1)) *
639  FEHermite<1>::hermite_raw_shape_deriv(bases1D[2],p(2));
640  break;
641  case 5:
642  return coef *
643  FEHermite<1>::hermite_raw_shape(bases1D[0],p(0)) *
644  FEHermite<1>::hermite_raw_shape(bases1D[1],p(1)) *
646  break;
647  default:
648  libmesh_error_msg("Invalid shape function derivative j = " << j);
649  }
650 
651  }
652  default:
653  libmesh_error_msg("ERROR: Unsupported element type " << type);
654  }
655  }
656  // by default throw an error
657  default:
658  libmesh_error_msg("ERROR: Unsupported polynomial order " << totalorder);
659  }
660 }
661 
662 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
663 
664 } // namespace libMesh
libMesh::FE::shape_second_deriv
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::FEInterface::n_shape_functions
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:446
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::FEInterface::shape_deriv_function
static shape_deriv_ptr shape_deriv_function(const unsigned int dim, const FEType &fe_t)
Definition: fe_interface.C:916
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::FEHermite::hermite_raw_shape_second_deriv
static Real hermite_raw_shape_second_deriv(const unsigned int basis_num, const Real xi)
1D hermite functions on unit interval
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::square_number_column
const unsigned char square_number_column[]
Definition: number_lookups.C:56
libMesh::cube_number_page
const unsigned char cube_number_page[]
Definition: number_lookups.C:308
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::cube_number_column
const unsigned char cube_number_column[]
Definition: number_lookups.C:84
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::cube_number_row
const unsigned char cube_number_row[]
Definition: number_lookups.C:196
libMesh::FE::shape_deriv
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::FEHermite::hermite_raw_shape
static Real hermite_raw_shape(const unsigned int basis_num, const Real xi)
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::Elem::default_order
virtual Order default_order() const =0
libMesh::FEMap::map_fe_type
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:47
libMesh::FEInterface::shape_deriv_ptr
Real(* shape_deriv_ptr)(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
Definition: fe_interface.h:345
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::FEFamily
FEFamily
Definition: enum_fe_family.h:34
libMesh::square_number_row
const unsigned char square_number_row[]
Definition: number_lookups.C:69
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::THIRD
Definition: enum_order.h:44
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::FEHermite::hermite_raw_shape_deriv
static Real hermite_raw_shape_deriv(const unsigned int basis_num, const Real xi)
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33