libMesh
fe_abstract.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 
20 // Local includes
21 #include "libmesh/fe.h"
22 #include "libmesh/libmesh_logging.h"
23 #include "libmesh/enum_elem_type.h"
24 
25 // For projection code:
26 #include "libmesh/boundary_info.h"
27 #include "libmesh/mesh_base.h"
28 #include "libmesh/dense_matrix.h"
29 #include "libmesh/dense_vector.h"
30 #include "libmesh/dof_map.h"
31 #include "libmesh/elem.h"
32 #include "libmesh/fe_interface.h"
33 #include "libmesh/numeric_vector.h"
34 #include "libmesh/periodic_boundaries.h"
35 #include "libmesh/periodic_boundary.h"
36 #include "libmesh/quadrature.h"
37 #include "libmesh/quadrature_gauss.h"
38 #include "libmesh/remote_elem.h"
39 #include "libmesh/tensor_value.h"
40 #include "libmesh/threads.h"
41 #include "libmesh/enum_elem_type.h"
42 
43 namespace libMesh
44 {
45 
46 FEAbstract::FEAbstract(const unsigned int d,
47  const FEType & fet) :
48  _fe_map( FEMap::build(fet) ),
49  dim(d),
50  calculations_started(false),
51  calculate_map(false),
52  calculate_phi(false),
53  calculate_dphi(false),
54  calculate_d2phi(false),
55  calculate_curl_phi(false),
56  calculate_div_phi(false),
57  calculate_dphiref(false),
58  fe_type(fet),
59  elem_type(INVALID_ELEM),
60  _p_level(0),
61  qrule(nullptr),
62  shapes_on_quadrature(false)
63 {
64 }
65 
66 
68 {
69 }
70 
71 
72 std::unique_ptr<FEAbstract> FEAbstract::build(const unsigned int dim,
73  const FEType & fet)
74 {
75  switch (dim)
76  {
77  // 0D
78  case 0:
79  {
80  switch (fet.family)
81  {
82  case CLOUGH:
83  return libmesh_make_unique<FE<0,CLOUGH>>(fet);
84 
85  case HERMITE:
86  return libmesh_make_unique<FE<0,HERMITE>>(fet);
87 
88  case LAGRANGE:
89  return libmesh_make_unique<FE<0,LAGRANGE>>(fet);
90 
91  case LAGRANGE_VEC:
92  return libmesh_make_unique<FE<0,LAGRANGE_VEC>>(fet);
93 
94  case L2_LAGRANGE:
95  return libmesh_make_unique<FE<0,L2_LAGRANGE>>(fet);
96 
97  case HIERARCHIC:
98  return libmesh_make_unique<FE<0,HIERARCHIC>>(fet);
99 
100  case L2_HIERARCHIC:
101  return libmesh_make_unique<FE<0,L2_HIERARCHIC>>(fet);
102 
103  case MONOMIAL:
104  return libmesh_make_unique<FE<0,MONOMIAL>>(fet);
105 
106  case MONOMIAL_VEC:
107  return libmesh_make_unique<FE<0,MONOMIAL_VEC>>(fet);
108 
109 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
110  case SZABAB:
111  return libmesh_make_unique<FE<0,SZABAB>>(fet);
112 
113  case BERNSTEIN:
114  return libmesh_make_unique<FE<0,BERNSTEIN>>(fet);
115 
116  case RATIONAL_BERNSTEIN:
117  return libmesh_make_unique<FE<0,RATIONAL_BERNSTEIN>>(fet);
118 #endif
119 
120  case XYZ:
121  return libmesh_make_unique<FEXYZ<0>>(fet);
122 
123  case SCALAR:
124  return libmesh_make_unique<FEScalar<0>>(fet);
125 
126  default:
127  libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
128  }
129  }
130  // 1D
131  case 1:
132  {
133  switch (fet.family)
134  {
135  case CLOUGH:
136  return libmesh_make_unique<FE<1,CLOUGH>>(fet);
137 
138  case HERMITE:
139  return libmesh_make_unique<FE<1,HERMITE>>(fet);
140 
141  case LAGRANGE:
142  return libmesh_make_unique<FE<1,LAGRANGE>>(fet);
143 
144  case LAGRANGE_VEC:
145  return libmesh_make_unique<FE<1,LAGRANGE_VEC>>(fet);
146 
147  case L2_LAGRANGE:
148  return libmesh_make_unique<FE<1,L2_LAGRANGE>>(fet);
149 
150  case HIERARCHIC:
151  return libmesh_make_unique<FE<1,HIERARCHIC>>(fet);
152 
153  case L2_HIERARCHIC:
154  return libmesh_make_unique<FE<1,L2_HIERARCHIC>>(fet);
155 
156  case MONOMIAL:
157  return libmesh_make_unique<FE<1,MONOMIAL>>(fet);
158 
159  case MONOMIAL_VEC:
160  return libmesh_make_unique<FE<1,MONOMIAL_VEC>>(fet);
161 
162 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
163  case SZABAB:
164  return libmesh_make_unique<FE<1,SZABAB>>(fet);
165 
166  case BERNSTEIN:
167  return libmesh_make_unique<FE<1,BERNSTEIN>>(fet);
168 
169  case RATIONAL_BERNSTEIN:
170  return libmesh_make_unique<FE<1,RATIONAL_BERNSTEIN>>(fet);
171 #endif
172 
173  case XYZ:
174  return libmesh_make_unique<FEXYZ<1>>(fet);
175 
176  case SCALAR:
177  return libmesh_make_unique<FEScalar<1>>(fet);
178 
179  default:
180  libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
181  }
182  }
183 
184 
185  // 2D
186  case 2:
187  {
188  switch (fet.family)
189  {
190  case CLOUGH:
191  return libmesh_make_unique<FE<2,CLOUGH>>(fet);
192 
193  case HERMITE:
194  return libmesh_make_unique<FE<2,HERMITE>>(fet);
195 
196  case LAGRANGE:
197  return libmesh_make_unique<FE<2,LAGRANGE>>(fet);
198 
199  case LAGRANGE_VEC:
200  return libmesh_make_unique<FE<2,LAGRANGE_VEC>>(fet);
201 
202  case L2_LAGRANGE:
203  return libmesh_make_unique<FE<2,L2_LAGRANGE>>(fet);
204 
205  case HIERARCHIC:
206  return libmesh_make_unique<FE<2,HIERARCHIC>>(fet);
207 
208  case L2_HIERARCHIC:
209  return libmesh_make_unique<FE<2,L2_HIERARCHIC>>(fet);
210 
211  case MONOMIAL:
212  return libmesh_make_unique<FE<2,MONOMIAL>>(fet);
213 
214  case MONOMIAL_VEC:
215  return libmesh_make_unique<FE<2,MONOMIAL_VEC>>(fet);
216 
217 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
218  case SZABAB:
219  return libmesh_make_unique<FE<2,SZABAB>>(fet);
220 
221  case BERNSTEIN:
222  return libmesh_make_unique<FE<2,BERNSTEIN>>(fet);
223 
224  case RATIONAL_BERNSTEIN:
225  return libmesh_make_unique<FE<2,RATIONAL_BERNSTEIN>>(fet);
226 #endif
227 
228  case XYZ:
229  return libmesh_make_unique<FEXYZ<2>>(fet);
230 
231  case SCALAR:
232  return libmesh_make_unique<FEScalar<2>>(fet);
233 
234  case NEDELEC_ONE:
235  return libmesh_make_unique<FENedelecOne<2>>(fet);
236 
237  case SUBDIVISION:
238  return libmesh_make_unique<FESubdivision>(fet);
239 
240  default:
241  libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
242  }
243  }
244 
245 
246  // 3D
247  case 3:
248  {
249  switch (fet.family)
250  {
251  case CLOUGH:
252  libmesh_error_msg("ERROR: Clough-Tocher elements currently only support 1D and 2D");
253 
254  case HERMITE:
255  return libmesh_make_unique<FE<3,HERMITE>>(fet);
256 
257  case LAGRANGE:
258  return libmesh_make_unique<FE<3,LAGRANGE>>(fet);
259 
260  case LAGRANGE_VEC:
261  return libmesh_make_unique<FE<3,LAGRANGE_VEC>>(fet);
262 
263  case L2_LAGRANGE:
264  return libmesh_make_unique<FE<3,L2_LAGRANGE>>(fet);
265 
266  case HIERARCHIC:
267  return libmesh_make_unique<FE<3,HIERARCHIC>>(fet);
268 
269  case L2_HIERARCHIC:
270  return libmesh_make_unique<FE<3,L2_HIERARCHIC>>(fet);
271 
272  case MONOMIAL:
273  return libmesh_make_unique<FE<3,MONOMIAL>>(fet);
274 
275  case MONOMIAL_VEC:
276  return libmesh_make_unique<FE<3,MONOMIAL_VEC>>(fet);
277 
278 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES
279  case SZABAB:
280  return libmesh_make_unique<FE<3,SZABAB>>(fet);
281 
282  case BERNSTEIN:
283  return libmesh_make_unique<FE<3,BERNSTEIN>>(fet);
284 
285  case RATIONAL_BERNSTEIN:
286  return libmesh_make_unique<FE<3,RATIONAL_BERNSTEIN>>(fet);
287 #endif
288 
289  case XYZ:
290  return libmesh_make_unique<FEXYZ<3>>(fet);
291 
292  case SCALAR:
293  return libmesh_make_unique<FEScalar<3>>(fet);
294 
295  case NEDELEC_ONE:
296  return libmesh_make_unique<FENedelecOne<3>>(fet);
297 
298  default:
299  libmesh_error_msg("ERROR: Bad FEType.family= " << fet.family);
300  }
301  }
302 
303  default:
304  libmesh_error_msg("Invalid dimension dim = " << dim);
305  }
306 }
307 
308 void FEAbstract::get_refspace_nodes(const ElemType itemType, std::vector<Point> & nodes)
309 {
310  switch(itemType)
311  {
312  case EDGE2:
313  {
314  nodes.resize(2);
315  nodes[0] = Point (-1.,0.,0.);
316  nodes[1] = Point (1.,0.,0.);
317  return;
318  }
319  case EDGE3:
320  {
321  nodes.resize(3);
322  nodes[0] = Point (-1.,0.,0.);
323  nodes[1] = Point (1.,0.,0.);
324  nodes[2] = Point (0.,0.,0.);
325  return;
326  }
327  case TRI3:
328  case TRISHELL3:
329  {
330  nodes.resize(3);
331  nodes[0] = Point (0.,0.,0.);
332  nodes[1] = Point (1.,0.,0.);
333  nodes[2] = Point (0.,1.,0.);
334  return;
335  }
336  case TRI6:
337  {
338  nodes.resize(6);
339  nodes[0] = Point (0.,0.,0.);
340  nodes[1] = Point (1.,0.,0.);
341  nodes[2] = Point (0.,1.,0.);
342  nodes[3] = Point (.5,0.,0.);
343  nodes[4] = Point (.5,.5,0.);
344  nodes[5] = Point (0.,.5,0.);
345  return;
346  }
347  case QUAD4:
348  case QUADSHELL4:
349  {
350  nodes.resize(4);
351  nodes[0] = Point (-1.,-1.,0.);
352  nodes[1] = Point (1.,-1.,0.);
353  nodes[2] = Point (1.,1.,0.);
354  nodes[3] = Point (-1.,1.,0.);
355  return;
356  }
357  case QUAD8:
358  case QUADSHELL8:
359  {
360  nodes.resize(8);
361  nodes[0] = Point (-1.,-1.,0.);
362  nodes[1] = Point (1.,-1.,0.);
363  nodes[2] = Point (1.,1.,0.);
364  nodes[3] = Point (-1.,1.,0.);
365  nodes[4] = Point (0.,-1.,0.);
366  nodes[5] = Point (1.,0.,0.);
367  nodes[6] = Point (0.,1.,0.);
368  nodes[7] = Point (-1.,0.,0.);
369  return;
370  }
371  case QUAD9:
372  {
373  nodes.resize(9);
374  nodes[0] = Point (-1.,-1.,0.);
375  nodes[1] = Point (1.,-1.,0.);
376  nodes[2] = Point (1.,1.,0.);
377  nodes[3] = Point (-1.,1.,0.);
378  nodes[4] = Point (0.,-1.,0.);
379  nodes[5] = Point (1.,0.,0.);
380  nodes[6] = Point (0.,1.,0.);
381  nodes[7] = Point (-1.,0.,0.);
382  nodes[8] = Point (0.,0.,0.);
383  return;
384  }
385  case TET4:
386  {
387  nodes.resize(4);
388  nodes[0] = Point (0.,0.,0.);
389  nodes[1] = Point (1.,0.,0.);
390  nodes[2] = Point (0.,1.,0.);
391  nodes[3] = Point (0.,0.,1.);
392  return;
393  }
394  case TET10:
395  {
396  nodes.resize(10);
397  nodes[0] = Point (0.,0.,0.);
398  nodes[1] = Point (1.,0.,0.);
399  nodes[2] = Point (0.,1.,0.);
400  nodes[3] = Point (0.,0.,1.);
401  nodes[4] = Point (.5,0.,0.);
402  nodes[5] = Point (.5,.5,0.);
403  nodes[6] = Point (0.,.5,0.);
404  nodes[7] = Point (0.,0.,.5);
405  nodes[8] = Point (.5,0.,.5);
406  nodes[9] = Point (0.,.5,.5);
407  return;
408  }
409  case HEX8:
410  {
411  nodes.resize(8);
412  nodes[0] = Point (-1.,-1.,-1.);
413  nodes[1] = Point (1.,-1.,-1.);
414  nodes[2] = Point (1.,1.,-1.);
415  nodes[3] = Point (-1.,1.,-1.);
416  nodes[4] = Point (-1.,-1.,1.);
417  nodes[5] = Point (1.,-1.,1.);
418  nodes[6] = Point (1.,1.,1.);
419  nodes[7] = Point (-1.,1.,1.);
420  return;
421  }
422  case HEX20:
423  {
424  nodes.resize(20);
425  nodes[0] = Point (-1.,-1.,-1.);
426  nodes[1] = Point (1.,-1.,-1.);
427  nodes[2] = Point (1.,1.,-1.);
428  nodes[3] = Point (-1.,1.,-1.);
429  nodes[4] = Point (-1.,-1.,1.);
430  nodes[5] = Point (1.,-1.,1.);
431  nodes[6] = Point (1.,1.,1.);
432  nodes[7] = Point (-1.,1.,1.);
433  nodes[8] = Point (0.,-1.,-1.);
434  nodes[9] = Point (1.,0.,-1.);
435  nodes[10] = Point (0.,1.,-1.);
436  nodes[11] = Point (-1.,0.,-1.);
437  nodes[12] = Point (-1.,-1.,0.);
438  nodes[13] = Point (1.,-1.,0.);
439  nodes[14] = Point (1.,1.,0.);
440  nodes[15] = Point (-1.,1.,0.);
441  nodes[16] = Point (0.,-1.,1.);
442  nodes[17] = Point (1.,0.,1.);
443  nodes[18] = Point (0.,1.,1.);
444  nodes[19] = Point (-1.,0.,1.);
445  return;
446  }
447  case HEX27:
448  {
449  nodes.resize(27);
450  nodes[0] = Point (-1.,-1.,-1.);
451  nodes[1] = Point (1.,-1.,-1.);
452  nodes[2] = Point (1.,1.,-1.);
453  nodes[3] = Point (-1.,1.,-1.);
454  nodes[4] = Point (-1.,-1.,1.);
455  nodes[5] = Point (1.,-1.,1.);
456  nodes[6] = Point (1.,1.,1.);
457  nodes[7] = Point (-1.,1.,1.);
458  nodes[8] = Point (0.,-1.,-1.);
459  nodes[9] = Point (1.,0.,-1.);
460  nodes[10] = Point (0.,1.,-1.);
461  nodes[11] = Point (-1.,0.,-1.);
462  nodes[12] = Point (-1.,-1.,0.);
463  nodes[13] = Point (1.,-1.,0.);
464  nodes[14] = Point (1.,1.,0.);
465  nodes[15] = Point (-1.,1.,0.);
466  nodes[16] = Point (0.,-1.,1.);
467  nodes[17] = Point (1.,0.,1.);
468  nodes[18] = Point (0.,1.,1.);
469  nodes[19] = Point (-1.,0.,1.);
470  nodes[20] = Point (0.,0.,-1.);
471  nodes[21] = Point (0.,-1.,0.);
472  nodes[22] = Point (1.,0.,0.);
473  nodes[23] = Point (0.,1.,0.);
474  nodes[24] = Point (-1.,0.,0.);
475  nodes[25] = Point (0.,0.,1.);
476  nodes[26] = Point (0.,0.,0.);
477  return;
478  }
479  case PRISM6:
480  {
481  nodes.resize(6);
482  nodes[0] = Point (0.,0.,-1.);
483  nodes[1] = Point (1.,0.,-1.);
484  nodes[2] = Point (0.,1.,-1.);
485  nodes[3] = Point (0.,0.,1.);
486  nodes[4] = Point (1.,0.,1.);
487  nodes[5] = Point (0.,1.,1.);
488  return;
489  }
490  case PRISM15:
491  {
492  nodes.resize(15);
493  nodes[0] = Point (0.,0.,-1.);
494  nodes[1] = Point (1.,0.,-1.);
495  nodes[2] = Point (0.,1.,-1.);
496  nodes[3] = Point (0.,0.,1.);
497  nodes[4] = Point (1.,0.,1.);
498  nodes[5] = Point (0.,1.,1.);
499  nodes[6] = Point (.5,0.,-1.);
500  nodes[7] = Point (.5,.5,-1.);
501  nodes[8] = Point (0.,.5,-1.);
502  nodes[9] = Point (0.,0.,0.);
503  nodes[10] = Point (1.,0.,0.);
504  nodes[11] = Point (0.,1.,0.);
505  nodes[12] = Point (.5,0.,1.);
506  nodes[13] = Point (.5,.5,1.);
507  nodes[14] = Point (0.,.5,1.);
508  return;
509  }
510  case PRISM18:
511  {
512  nodes.resize(18);
513  nodes[0] = Point (0.,0.,-1.);
514  nodes[1] = Point (1.,0.,-1.);
515  nodes[2] = Point (0.,1.,-1.);
516  nodes[3] = Point (0.,0.,1.);
517  nodes[4] = Point (1.,0.,1.);
518  nodes[5] = Point (0.,1.,1.);
519  nodes[6] = Point (.5,0.,-1.);
520  nodes[7] = Point (.5,.5,-1.);
521  nodes[8] = Point (0.,.5,-1.);
522  nodes[9] = Point (0.,0.,0.);
523  nodes[10] = Point (1.,0.,0.);
524  nodes[11] = Point (0.,1.,0.);
525  nodes[12] = Point (.5,0.,1.);
526  nodes[13] = Point (.5,.5,1.);
527  nodes[14] = Point (0.,.5,1.);
528  nodes[15] = Point (.5,0.,0.);
529  nodes[16] = Point (.5,.5,0.);
530  nodes[17] = Point (0.,.5,0.);
531  return;
532  }
533  case PYRAMID5:
534  {
535  nodes.resize(5);
536  nodes[0] = Point (-1.,-1.,0.);
537  nodes[1] = Point (1.,-1.,0.);
538  nodes[2] = Point (1.,1.,0.);
539  nodes[3] = Point (-1.,1.,0.);
540  nodes[4] = Point (0.,0.,1.);
541  return;
542  }
543  case PYRAMID13:
544  {
545  nodes.resize(13);
546 
547  // base corners
548  nodes[0] = Point (-1.,-1.,0.);
549  nodes[1] = Point (1.,-1.,0.);
550  nodes[2] = Point (1.,1.,0.);
551  nodes[3] = Point (-1.,1.,0.);
552 
553  // apex
554  nodes[4] = Point (0.,0.,1.);
555 
556  // base midedge
557  nodes[5] = Point (0.,-1.,0.);
558  nodes[6] = Point (1.,0.,0.);
559  nodes[7] = Point (0.,1.,0.);
560  nodes[8] = Point (-1,0.,0.);
561 
562  // lateral midedge
563  nodes[9] = Point (-.5,-.5,.5);
564  nodes[10] = Point (.5,-.5,.5);
565  nodes[11] = Point (.5,.5,.5);
566  nodes[12] = Point (-.5,.5,.5);
567 
568  return;
569  }
570  case PYRAMID14:
571  {
572  nodes.resize(14);
573 
574  // base corners
575  nodes[0] = Point (-1.,-1.,0.);
576  nodes[1] = Point (1.,-1.,0.);
577  nodes[2] = Point (1.,1.,0.);
578  nodes[3] = Point (-1.,1.,0.);
579 
580  // apex
581  nodes[4] = Point (0.,0.,1.);
582 
583  // base midedge
584  nodes[5] = Point (0.,-1.,0.);
585  nodes[6] = Point (1.,0.,0.);
586  nodes[7] = Point (0.,1.,0.);
587  nodes[8] = Point (-1,0.,0.);
588 
589  // lateral midedge
590  nodes[9] = Point (-.5,-.5,.5);
591  nodes[10] = Point (.5,-.5,.5);
592  nodes[11] = Point (.5,.5,.5);
593  nodes[12] = Point (-.5,.5,.5);
594 
595  // base center
596  nodes[13] = Point (0.,0.,0.);
597 
598  return;
599  }
600 
601  default:
602  libmesh_error_msg("ERROR: Unknown element type " << itemType);
603  }
604 }
605 
606 bool FEAbstract::on_reference_element(const Point & p, const ElemType t, const Real eps)
607 {
608  libmesh_assert_greater_equal (eps, 0.);
609 
610  const Real xi = p(0);
611 #if LIBMESH_DIM > 1
612  const Real eta = p(1);
613 #else
614  const Real eta = 0.;
615 #endif
616 #if LIBMESH_DIM > 2
617  const Real zeta = p(2);
618 #else
619  const Real zeta = 0.;
620 #endif
621 
622  switch (t)
623  {
624  case NODEELEM:
625  {
626  return (!xi && !eta && !zeta);
627  }
628  case EDGE2:
629  case EDGE3:
630  case EDGE4:
631  {
632  // The reference 1D element is [-1,1].
633  if ((xi >= -1.-eps) &&
634  (xi <= 1.+eps))
635  return true;
636 
637  return false;
638  }
639 
640 
641  case TRI3:
642  case TRISHELL3:
643  case TRI6:
644  {
645  // The reference triangle is isosceles
646  // and is bound by xi=0, eta=0, and xi+eta=1.
647  if ((xi >= 0.-eps) &&
648  (eta >= 0.-eps) &&
649  ((xi + eta) <= 1.+eps))
650  return true;
651 
652  return false;
653  }
654 
655 
656  case QUAD4:
657  case QUADSHELL4:
658  case QUAD8:
659  case QUADSHELL8:
660  case QUAD9:
661  {
662  // The reference quadrilateral element is [-1,1]^2.
663  if ((xi >= -1.-eps) &&
664  (xi <= 1.+eps) &&
665  (eta >= -1.-eps) &&
666  (eta <= 1.+eps))
667  return true;
668 
669  return false;
670  }
671 
672 
673  case TET4:
674  case TET10:
675  {
676  // The reference tetrahedral is isosceles
677  // and is bound by xi=0, eta=0, zeta=0,
678  // and xi+eta+zeta=1.
679  if ((xi >= 0.-eps) &&
680  (eta >= 0.-eps) &&
681  (zeta >= 0.-eps) &&
682  ((xi + eta + zeta) <= 1.+eps))
683  return true;
684 
685  return false;
686  }
687 
688 
689  case HEX8:
690  case HEX20:
691  case HEX27:
692  {
693  /*
694  if ((xi >= -1.) &&
695  (xi <= 1.) &&
696  (eta >= -1.) &&
697  (eta <= 1.) &&
698  (zeta >= -1.) &&
699  (zeta <= 1.))
700  return true;
701  */
702 
703  // The reference hexahedral element is [-1,1]^3.
704  if ((xi >= -1.-eps) &&
705  (xi <= 1.+eps) &&
706  (eta >= -1.-eps) &&
707  (eta <= 1.+eps) &&
708  (zeta >= -1.-eps) &&
709  (zeta <= 1.+eps))
710  {
711  // libMesh::out << "Strange Point:\n";
712  // p.print();
713  return true;
714  }
715 
716  return false;
717  }
718 
719  case PRISM6:
720  case PRISM15:
721  case PRISM18:
722  {
723  // Figure this one out...
724  // inside the reference triangle with zeta in [-1,1]
725  if ((xi >= 0.-eps) &&
726  (eta >= 0.-eps) &&
727  (zeta >= -1.-eps) &&
728  (zeta <= 1.+eps) &&
729  ((xi + eta) <= 1.+eps))
730  return true;
731 
732  return false;
733  }
734 
735 
736  case PYRAMID5:
737  case PYRAMID13:
738  case PYRAMID14:
739  {
740  // Check that the point is on the same side of all the faces
741  // by testing whether:
742  //
743  // n_i.(x - x_i) <= 0
744  //
745  // for each i, where:
746  // n_i is the outward normal of face i,
747  // x_i is a point on face i.
748  if ((-eta - 1. + zeta <= 0.+eps) &&
749  ( xi - 1. + zeta <= 0.+eps) &&
750  ( eta - 1. + zeta <= 0.+eps) &&
751  ( -xi - 1. + zeta <= 0.+eps) &&
752  ( zeta >= 0.-eps))
753  return true;
754 
755  return false;
756  }
757 
758 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
759  case INFHEX8:
760  case INFHEX16:
761  case INFHEX18:
762  {
763  // The reference infhex8 is a [-1,1]^3.
764  if ((xi >= -1.-eps) &&
765  (xi <= 1.+eps) &&
766  (eta >= -1.-eps) &&
767  (eta <= 1.+eps) &&
768  (zeta >= -1.-eps) &&
769  (zeta <= 1.+eps))
770  {
771  return true;
772  }
773  return false;
774  }
775 
776  case INFPRISM6:
777  case INFPRISM12:
778  {
779  // inside the reference triangle with zeta in [-1,1]
780  if ((xi >= 0.-eps) &&
781  (eta >= 0.-eps) &&
782  (zeta >= -1.-eps) &&
783  (zeta <= 1.+eps) &&
784  ((xi + eta) <= 1.+eps))
785  {
786  return true;
787  }
788 
789  return false;
790  }
791 #endif
792 
793  default:
794  libmesh_error_msg("ERROR: Unknown element type " << t);
795  }
796 
797  // If we get here then the point is _not_ in the
798  // reference element. Better return false.
799 
800  return false;
801 }
802 
803 
804 
805 void FEAbstract::print_JxW(std::ostream & os) const
806 {
807  this->_fe_map->print_JxW(os);
808 }
809 
810 
811 
812 void FEAbstract::print_xyz(std::ostream & os) const
813 {
814  this->_fe_map->print_xyz(os);
815 }
816 
817 
818 void FEAbstract::print_info(std::ostream & os) const
819 {
820  os << "phi[i][j]: Shape function i at quadrature pt. j" << std::endl;
821  this->print_phi(os);
822 
823  os << "dphi[i][j]: Shape function i's gradient at quadrature pt. j" << std::endl;
824  this->print_dphi(os);
825 
826  os << "XYZ locations of the quadrature pts." << std::endl;
827  this->print_xyz(os);
828 
829  os << "Values of JxW at the quadrature pts." << std::endl;
830  this->print_JxW(os);
831 }
832 
833 
834 std::ostream & operator << (std::ostream & os, const FEAbstract & fe)
835 {
836  fe.print_info(os);
837  return os;
838 }
839 
840 
841 
842 #ifdef LIBMESH_ENABLE_AMR
843 
844 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
846  const Elem * elem)
847 {
848  libmesh_assert(elem);
849 
850  const unsigned int Dim = elem->dim();
851 
852  // Only constrain elements in 2,3D.
853  if (Dim == 1)
854  return;
855 
856  // Only constrain active and ancestor elements
857  if (elem->subactive())
858  return;
859 
860  const FEFamily mapping_family = FEMap::map_fe_type(*elem);
861  const FEType fe_type(elem->default_order(), mapping_family);
862 
863  // Pull objects out of the loop to reduce heap operations
864  std::vector<const Node *> my_nodes, parent_nodes;
865  std::unique_ptr<const Elem> my_side, parent_side;
866 
867  // Look at the element faces. Check to see if we need to
868  // build constraints.
869  for (auto s : elem->side_index_range())
870  if (elem->neighbor_ptr(s) != nullptr &&
871  elem->neighbor_ptr(s) != remote_elem)
872  if (elem->neighbor_ptr(s)->level() < elem->level()) // constrain dofs shared between
873  { // this element and ones coarser
874  // than this element.
875  // Get pointers to the elements of interest and its parent.
876  const Elem * parent = elem->parent();
877 
878  // This can't happen... Only level-0 elements have nullptr
879  // parents, and no level-0 elements can be at a higher
880  // level than their neighbors!
881  libmesh_assert(parent);
882 
883  elem->build_side_ptr(my_side, s);
884  parent->build_side_ptr(parent_side, s);
885 
886  const unsigned int n_side_nodes = my_side->n_nodes();
887 
888  my_nodes.clear();
889  my_nodes.reserve (n_side_nodes);
890  parent_nodes.clear();
891  parent_nodes.reserve (n_side_nodes);
892 
893  for (unsigned int n=0; n != n_side_nodes; ++n)
894  my_nodes.push_back(my_side->node_ptr(n));
895 
896  for (unsigned int n=0; n != n_side_nodes; ++n)
897  parent_nodes.push_back(parent_side->node_ptr(n));
898 
899  for (unsigned int my_side_n=0;
900  my_side_n < n_side_nodes;
901  my_side_n++)
902  {
903  libmesh_assert_less (my_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));
904 
905  const Node * my_node = my_nodes[my_side_n];
906 
907  // The support point of the DOF
908  const Point & support_point = *my_node;
909 
910  // Figure out where my node lies on their reference element.
911  const Point mapped_point = FEMap::inverse_map(Dim-1,
912  parent_side.get(),
913  support_point);
914 
915  // Compute the parent's side shape function values.
916  for (unsigned int their_side_n=0;
917  their_side_n < n_side_nodes;
918  their_side_n++)
919  {
920  libmesh_assert_less (their_side_n, FEInterface::n_dofs(Dim-1, fe_type, parent_side->type()));
921 
922  const Node * their_node = parent_nodes[their_side_n];
923  libmesh_assert(their_node);
924 
925  const Real their_value = FEInterface::shape(Dim-1,
926  fe_type,
927  parent_side->type(),
928  their_side_n,
929  mapped_point);
930 
931  const Real their_mag = std::abs(their_value);
932 #ifdef DEBUG
933  // Protect for the case u_i ~= u_j,
934  // in which case i better equal j.
935  if (their_mag > 0.999)
936  {
937  libmesh_assert_equal_to (my_node, their_node);
938  libmesh_assert_less (std::abs(their_value - 1.), 0.001);
939  }
940  else
941 #endif
942  // To make nodal constraints useful for constructing
943  // sparsity patterns faster, we need to get EVERY
944  // POSSIBLE constraint coupling identified, even if
945  // there is no coupling in the isoparametric
946  // Lagrange case.
947  if (their_mag < 1.e-5)
948  {
949  // since we may be running this method concurrently
950  // on multiple threads we need to acquire a lock
951  // before modifying the shared constraint_row object.
952  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
953 
954  // A reference to the constraint row.
955  NodeConstraintRow & constraint_row = constraints[my_node].first;
956 
957  constraint_row.insert(std::make_pair (their_node,
958  0.));
959  }
960  // To get nodal coordinate constraints right, only
961  // add non-zero and non-identity values for Lagrange
962  // basis functions.
963  else // (1.e-5 <= their_mag <= .999)
964  {
965  // since we may be running this method concurrently
966  // on multiple threads we need to acquire a lock
967  // before modifying the shared constraint_row object.
968  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
969 
970  // A reference to the constraint row.
971  NodeConstraintRow & constraint_row = constraints[my_node].first;
972 
973  constraint_row.insert(std::make_pair (their_node,
974  their_value));
975  }
976  }
977  }
978  }
979 }
980 
981 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
982 
983 #endif // #ifdef LIBMESH_ENABLE_AMR
984 
985 
986 
987 #ifdef LIBMESH_ENABLE_PERIODIC
988 
989 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS
991  const PeriodicBoundaries & boundaries,
992  const MeshBase & mesh,
993  const PointLocatorBase * point_locator,
994  const Elem * elem)
995 {
996  // Only bother if we truly have periodic boundaries
997  if (boundaries.empty())
998  return;
999 
1000  libmesh_assert(elem);
1001 
1002  // Only constrain active elements with this method
1003  if (!elem->active())
1004  return;
1005 
1006  const unsigned int Dim = elem->dim();
1007 
1008  const FEFamily mapping_family = FEMap::map_fe_type(*elem);
1009  const FEType fe_type(elem->default_order(), mapping_family);
1010 
1011  // Pull objects out of the loop to reduce heap operations
1012  std::vector<const Node *> my_nodes, neigh_nodes;
1013  std::unique_ptr<const Elem> my_side, neigh_side;
1014 
1015  // Look at the element faces. Check to see if we need to
1016  // build constraints.
1017  std::vector<boundary_id_type> bc_ids;
1018  for (auto s : elem->side_index_range())
1019  {
1020  if (elem->neighbor_ptr(s))
1021  continue;
1022 
1023  mesh.get_boundary_info().boundary_ids (elem, s, bc_ids);
1024  for (const auto & boundary_id : bc_ids)
1025  {
1026  const PeriodicBoundaryBase * periodic = boundaries.boundary(boundary_id);
1027  if (periodic)
1028  {
1029  libmesh_assert(point_locator);
1030 
1031  // Get pointers to the element's neighbor.
1032  const Elem * neigh = boundaries.neighbor(boundary_id, *point_locator, elem, s);
1033 
1034  // h refinement constraints:
1035  // constrain dofs shared between
1036  // this element and ones as coarse
1037  // as or coarser than this element.
1038  if (neigh->level() <= elem->level())
1039  {
1040  unsigned int s_neigh =
1041  mesh.get_boundary_info().side_with_boundary_id(neigh, periodic->pairedboundary);
1042  libmesh_assert_not_equal_to (s_neigh, libMesh::invalid_uint);
1043 
1044 #ifdef LIBMESH_ENABLE_AMR
1045  libmesh_assert(neigh->active());
1046 #endif // #ifdef LIBMESH_ENABLE_AMR
1047 
1048  elem->build_side_ptr(my_side, s);
1049  neigh->build_side_ptr(neigh_side, s_neigh);
1050 
1051  const unsigned int n_side_nodes = my_side->n_nodes();
1052 
1053  my_nodes.clear();
1054  my_nodes.reserve (n_side_nodes);
1055  neigh_nodes.clear();
1056  neigh_nodes.reserve (n_side_nodes);
1057 
1058  for (unsigned int n=0; n != n_side_nodes; ++n)
1059  my_nodes.push_back(my_side->node_ptr(n));
1060 
1061  for (unsigned int n=0; n != n_side_nodes; ++n)
1062  neigh_nodes.push_back(neigh_side->node_ptr(n));
1063 
1064  // Make sure we're not adding recursive constraints
1065  // due to the redundancy in the way we add periodic
1066  // boundary constraints, or adding constraints to
1067  // nodes that already have AMR constraints
1068  std::vector<bool> skip_constraint(n_side_nodes, false);
1069 
1070  for (unsigned int my_side_n=0;
1071  my_side_n < n_side_nodes;
1072  my_side_n++)
1073  {
1074  libmesh_assert_less (my_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));
1075 
1076  const Node * my_node = my_nodes[my_side_n];
1077 
1078  // Figure out where my node lies on their reference element.
1079  const Point neigh_point = periodic->get_corresponding_pos(*my_node);
1080 
1081  const Point mapped_point =
1082  FEMap::inverse_map(Dim-1, neigh_side.get(),
1083  neigh_point);
1084 
1085  // If we've already got a constraint on this
1086  // node, then the periodic constraint is
1087  // redundant
1088  {
1089  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
1090 
1091  if (constraints.count(my_node))
1092  {
1093  skip_constraint[my_side_n] = true;
1094  continue;
1095  }
1096  }
1097 
1098  // Compute the neighbors's side shape function values.
1099  for (unsigned int their_side_n=0;
1100  their_side_n < n_side_nodes;
1101  their_side_n++)
1102  {
1103  libmesh_assert_less (their_side_n, FEInterface::n_dofs(Dim-1, fe_type, neigh_side->type()));
1104 
1105  const Node * their_node = neigh_nodes[their_side_n];
1106 
1107  // If there's a constraint on an opposing node,
1108  // we need to see if it's constrained by
1109  // *our side* making any periodic constraint
1110  // on us recursive
1111  {
1112  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
1113 
1114  if (!constraints.count(their_node))
1115  continue;
1116 
1117  const NodeConstraintRow & their_constraint_row =
1118  constraints[their_node].first;
1119 
1120  for (unsigned int orig_side_n=0;
1121  orig_side_n < n_side_nodes;
1122  orig_side_n++)
1123  {
1124  libmesh_assert_less (orig_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));
1125 
1126  const Node * orig_node = my_nodes[orig_side_n];
1127 
1128  if (their_constraint_row.count(orig_node))
1129  skip_constraint[orig_side_n] = true;
1130  }
1131  }
1132  }
1133  }
1134  for (unsigned int my_side_n=0;
1135  my_side_n < n_side_nodes;
1136  my_side_n++)
1137  {
1138  libmesh_assert_less (my_side_n, FEInterface::n_dofs(Dim-1, fe_type, my_side->type()));
1139 
1140  if (skip_constraint[my_side_n])
1141  continue;
1142 
1143  const Node * my_node = my_nodes[my_side_n];
1144 
1145  // Figure out where my node lies on their reference element.
1146  const Point neigh_point = periodic->get_corresponding_pos(*my_node);
1147 
1148  // Figure out where my node lies on their reference element.
1149  const Point mapped_point =
1150  FEMap::inverse_map(Dim-1, neigh_side.get(),
1151  neigh_point);
1152 
1153  for (unsigned int their_side_n=0;
1154  their_side_n < n_side_nodes;
1155  their_side_n++)
1156  {
1157  libmesh_assert_less (their_side_n, FEInterface::n_dofs(Dim-1, fe_type, neigh_side->type()));
1158 
1159  const Node * their_node = neigh_nodes[their_side_n];
1160  libmesh_assert(their_node);
1161 
1162  const Real their_value = FEInterface::shape(Dim-1,
1163  fe_type,
1164  neigh_side->type(),
1165  their_side_n,
1166  mapped_point);
1167 
1168  // since we may be running this method concurrently
1169  // on multiple threads we need to acquire a lock
1170  // before modifying the shared constraint_row object.
1171  {
1172  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
1173 
1174  NodeConstraintRow & constraint_row =
1175  constraints[my_node].first;
1176 
1177  constraint_row.insert(std::make_pair(their_node,
1178  their_value));
1179  }
1180  }
1181  }
1182  }
1183  }
1184  }
1185  }
1186 }
1187 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS
1188 
1189 #endif // LIBMESH_ENABLE_PERIODIC
1190 
1191 
1192 } // namespace libMesh
libMesh::HEX20
Definition: enum_elem_type.h:48
libMesh::PRISM6
Definition: enum_elem_type.h:50
libMesh::PeriodicBoundaries
We're using a class instead of a typedef to allow forward declarations and future flexibility.
Definition: periodic_boundaries.h:51
libMesh::invalid_uint
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value.
Definition: libmesh.h:249
libMesh::PeriodicBoundaryBase::get_corresponding_pos
virtual Point get_corresponding_pos(const Point &pt) const =0
This function should be overridden by derived classes to define how one finds corresponding nodes on ...
libMesh::CLOUGH
Definition: enum_fe_family.h:53
libMesh::Elem::build_side_ptr
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=true)=0
libMesh::FEType::family
FEFamily family
The type of finite element.
Definition: fe_type.h:203
libMesh::Elem::level
unsigned int level() const
Definition: elem.h:2478
libMesh::L2_HIERARCHIC
Definition: enum_fe_family.h:40
libMesh::HEX8
Definition: enum_elem_type.h:47
libMesh::FEAbstract::print_xyz
void print_xyz(std::ostream &os) const
Prints the spatial location of each quadrature point (on the physical element).
Definition: fe_abstract.C:812
libMesh::EDGE4
Definition: enum_elem_type.h:37
libMesh::INFHEX8
Definition: enum_elem_type.h:60
libMesh::FEAbstract::fe_type
FEType fe_type
The finite element type for this object.
Definition: fe_abstract.h:585
libMesh::HERMITE
Definition: enum_fe_family.h:54
libMesh::PeriodicBoundaryBase
The base class for defining periodic boundaries.
Definition: periodic_boundary_base.h:48
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Elem::dim
virtual unsigned short dim() const =0
libMesh::operator<<
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
Definition: fe_type.h:164
libMesh::TET10
Definition: enum_elem_type.h:46
libMesh::XYZ
Definition: enum_fe_family.h:46
libMesh::FEAbstract
This class forms the foundation from which generic finite elements may be derived.
Definition: fe_abstract.h:100
libMesh::FEAbstract::print_dphi
virtual void print_dphi(std::ostream &os) const =0
Prints the value of each shape function's derivative at each quadrature point.
libMesh::MONOMIAL_VEC
Definition: enum_fe_family.h:62
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::Elem::active
bool active() const
Definition: elem.h:2345
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::TET4
Definition: enum_elem_type.h:45
libMesh::PeriodicBoundaries::boundary
PeriodicBoundaryBase * boundary(boundary_id_type id)
Definition: periodic_boundaries.C:40
libMesh::PRISM15
Definition: enum_elem_type.h:51
libMesh::SZABAB
Definition: enum_fe_family.h:44
libMesh::FEAbstract::print_info
void print_info(std::ostream &os) const
Prints all the relevant information about the current element.
Definition: fe_abstract.C:818
libMesh::INFPRISM6
Definition: enum_elem_type.h:63
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::Elem::subactive
bool subactive() const
Definition: elem.h:2363
libMesh::MeshBase
This is the MeshBase class.
Definition: mesh_base.h:78
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::HEX27
Definition: enum_elem_type.h:49
libMesh::FEMap::inverse_map
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
Definition: fe_map.C:1622
libMesh::FEAbstract::print_phi
virtual void print_phi(std::ostream &os) const =0
Prints the value of each shape function at each quadrature point.
libMesh::Threads::spin_mtx
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:29
libMesh::FEAbstract::build
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
Definition: fe_abstract.C:72
libMesh::INFHEX18
Definition: enum_elem_type.h:62
libMesh::QUAD4
Definition: enum_elem_type.h:41
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::TRI3
Definition: enum_elem_type.h:39
libMesh::FEAbstract::_fe_map
std::unique_ptr< FEMap > _fe_map
Definition: fe_abstract.h:524
libMesh::Node
A Node is like a Point, but with more information.
Definition: node.h:52
libMesh::Elem::default_order
virtual Order default_order() const =0
libMesh::FEInterface::n_dofs
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:472
libMesh::FEMap::map_fe_type
static FEFamily map_fe_type(const Elem &elem)
Definition: fe_map.C:47
libMesh::RATIONAL_BERNSTEIN
Definition: enum_fe_family.h:64
libMesh::BERNSTEIN
Definition: enum_fe_family.h:43
libMesh::INFHEX16
Definition: enum_elem_type.h:61
libMesh::FEAbstract::dim
const unsigned int dim
The dimensionality of the object.
Definition: fe_abstract.h:530
libMesh::FEAbstract::compute_node_constraints
static void compute_node_constraints(NodeConstraints &constraints, const Elem *elem)
Computes the nodal constraint contributions (for non-conforming adapted meshes), using Lagrange geome...
Definition: fe_abstract.C:845
libMesh::FEAbstract::compute_periodic_node_constraints
static void compute_periodic_node_constraints(NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem)
Computes the node position constraint equation contributions (for meshes with periodic boundary condi...
Definition: fe_abstract.C:990
libMesh::HIERARCHIC
Definition: enum_fe_family.h:37
libMesh::MONOMIAL
Definition: enum_fe_family.h:39
libMesh::INFPRISM12
Definition: enum_elem_type.h:64
libMesh::Elem::parent
const Elem * parent() const
Definition: elem.h:2434
libMesh::FEAbstract::on_reference_element
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
Definition: fe_abstract.C:606
libMesh::TRI6
Definition: enum_elem_type.h:40
libMesh::INVALID_ELEM
Definition: enum_elem_type.h:75
libMesh::NodeConstraints
The Node constraint storage format.
Definition: dof_map.h:153
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::PeriodicBoundaries::neighbor
const Elem * neighbor(boundary_id_type boundary_id, const PointLocatorBase &point_locator, const Elem *e, unsigned int side) const
Definition: periodic_boundaries.C:61
libMesh::FEFamily
FEFamily
Definition: enum_fe_family.h:34
libMesh::QUADSHELL8
Definition: enum_elem_type.h:73
libMesh::PYRAMID5
Definition: enum_elem_type.h:53
libMesh::NEDELEC_ONE
Definition: enum_fe_family.h:61
libMesh::EDGE3
Definition: enum_elem_type.h:36
libMesh::QUADSHELL4
Definition: enum_elem_type.h:72
libMesh::FEAbstract::~FEAbstract
virtual ~FEAbstract()
Destructor.
Definition: fe_abstract.C:67
libMesh::Elem::side_index_range
IntRange< unsigned short > side_index_range() const
Definition: elem.h:2188
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::NODEELEM
Definition: enum_elem_type.h:66
libMesh::FEAbstract::get_refspace_nodes
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:308
libMesh::L2_LAGRANGE
Definition: enum_fe_family.h:41
libMesh::QUAD9
Definition: enum_elem_type.h:43
libMesh::FEMap
Class contained in FE that encapsulates mapping (i.e.
Definition: fe_map.h:47
libMesh::LAGRANGE_VEC
Definition: enum_fe_family.h:60
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::SCALAR
Definition: enum_fe_family.h:58
libMesh::SUBDIVISION
Definition: enum_fe_family.h:55
libMesh::Elem::neighbor_ptr
const Elem * neighbor_ptr(unsigned int i) const
Definition: elem.h:2085
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::PeriodicBoundaryBase::pairedboundary
boundary_id_type pairedboundary
Definition: periodic_boundary_base.h:58
libMesh::PRISM18
Definition: enum_elem_type.h:52
libMesh::TRISHELL3
Definition: enum_elem_type.h:71
libMesh::FEAbstract::FEAbstract
FEAbstract(const unsigned int dim, const FEType &fet)
Constructor.
Definition: fe_abstract.C:46
libMesh::PointLocatorBase
This is the base class for point locators.
Definition: point_locator_base.h:62
libMesh::remote_elem
const RemoteElem * remote_elem
Definition: remote_elem.C:57
libMesh::FEAbstract::print_JxW
void print_JxW(std::ostream &os) const
Prints the Jacobian times the weight for each quadrature point.
Definition: fe_abstract.C:805
libMesh::PYRAMID13
Definition: enum_elem_type.h:54
libMesh::PYRAMID14
Definition: enum_elem_type.h:55
libMesh::EDGE2
Definition: enum_elem_type.h:35
libMesh::QUAD8
Definition: enum_elem_type.h:42
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::NodeConstraintRow
std::map< const Node *, Real, std::less< const Node * >, Threads::scalable_allocator< std::pair< const Node *const, Real > > > NodeConstraintRow
A row of the Node constraint mapping.
Definition: dof_map.h:145
libMesh::FEInterface::shape
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
Definition: fe_interface.C:687