libMesh
inf_fe_static.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 // Local includes
20 #include "libmesh/libmesh_config.h"
21 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
22 #include "libmesh/inf_fe.h"
23 #include "libmesh/inf_fe_macro.h"
24 #include "libmesh/fe.h"
25 #include "libmesh/fe_interface.h"
26 #include "libmesh/fe_compute_data.h"
27 #include "libmesh/elem.h"
28 
29 namespace libMesh
30 {
31 
32 
33 // ------------------------------------------------------------
34 // InfFE class static member initialization
35 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
37 
38 #ifdef DEBUG
39 
40 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
42 
43 
44 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
46 
47 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
49 
50 #endif
51 
52 
53 
54 
55 // ------------------------------------------------------------
56 // InfFE static class members
57 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
58 unsigned int InfFE<Dim,T_radial,T_map>::n_dofs (const FEType & fet,
59  const ElemType inf_elem_type)
60 {
61  const ElemType base_et (InfFEBase::get_elem_type(inf_elem_type));
62 
63  if (Dim > 1)
64  return FEInterface::n_dofs(Dim-1, fet, base_et) *
66  else
68 }
69 
70 
71 
72 
73 
74 
75 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
77  const ElemType inf_elem_type,
78  const unsigned int n)
79 {
80  const ElemType base_et (InfFEBase::get_elem_type(inf_elem_type));
81 
82  unsigned int n_base, n_radial;
83  compute_node_indices(inf_elem_type, n, n_base, n_radial);
84 
85  // libMesh::out << "elem_type=" << inf_elem_type
86  // << ", fet.radial_order=" << fet.radial_order
87  // << ", n=" << n
88  // << ", n_radial=" << n_radial
89  // << ", n_base=" << n_base
90  // << std::endl;
91 
92  if (Dim > 1)
93  return FEInterface::n_dofs_at_node(Dim-1, fet, base_et, n_base)
95  else
96  return InfFERadial::n_dofs_at_node(fet.radial_order, n_radial);
97 }
98 
99 
100 
101 
102 
103 
104 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
106  const ElemType inf_elem_type)
107 {
108  const ElemType base_et (InfFEBase::get_elem_type(inf_elem_type));
109 
110  if (Dim > 1)
111  return FEInterface::n_dofs_per_elem(Dim-1, fet, base_et)
113  else
115 }
116 
117 
118 
119 
120 
121 
122 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
124  const Elem * /* elem */,
125  const std::vector<Number> & /* elem_soln */,
126  std::vector<Number> & nodal_soln)
127 {
128 #ifdef DEBUG
129  if (!_warned_for_nodal_soln)
130  {
131  libMesh::err << "WARNING: nodal_soln(...) does _not_ work for infinite elements." << std::endl
132  << " Will return an empty nodal solution. Use " << std::endl
133  << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!" << std::endl;
134  _warned_for_nodal_soln = true;
135  }
136 #endif
137 
138  /*
139  * In the base the infinite element couples to
140  * conventional finite elements. To not destroy
141  * the values there, clear \p nodal_soln. This
142  * indicates to the user of \p nodal_soln to
143  * not use this result.
144  */
145  nodal_soln.clear();
146  libmesh_assert (nodal_soln.empty());
147  return;
148 }
149 
150 
151 
152 
153 
154 
155 
156 
157 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
159  const ElemType inf_elem_type,
160  const unsigned int i,
161  const Point & p)
162 {
163  libmesh_assert_not_equal_to (Dim, 0);
164 
165 #ifdef DEBUG
166  // this makes only sense when used for mapping
167  if ((T_radial != INFINITE_MAP) && !_warned_for_shape)
168  {
169  libMesh::err << "WARNING: InfFE<Dim,T_radial,T_map>::shape(...) does _not_" << std::endl
170  << " return the correct trial function! Use " << std::endl
171  << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!"
172  << std::endl;
173  _warned_for_shape = true;
174  }
175 #endif
176 
177  const ElemType base_et (InfFEBase::get_elem_type(inf_elem_type));
178  const Order o_radial (fet.radial_order);
179  const Real v (p(Dim-1));
180 
181  unsigned int i_base, i_radial;
182  compute_shape_indices(fet, inf_elem_type, i, i_base, i_radial);
183 
184  //TODO:[SP/DD] exp(ikr) is still missing here!
185  // but is it intended? It would be probably somehow nice, but than it would be Number, not Real !
186  // --> thus it would destroy the interface...
187  if (Dim > 1)
188  return FEInterface::shape(Dim-1, fet, base_et, i_base, p)
189  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
190  * InfFERadial::decay(Dim,v);
191  else
192  return InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
193  * InfFERadial::decay(Dim,v);
194 }
195 
196 
197 
198 
199 
200 
201 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
203  const Elem * inf_elem,
204  const unsigned int i,
205  const Point & p)
206 {
207  libmesh_assert(inf_elem);
208  libmesh_assert_not_equal_to (Dim, 0);
209 
210 #ifdef DEBUG
211  // this makes only sense when used for mapping
212  if ((T_radial != INFINITE_MAP) && !_warned_for_shape)
213  {
214  libMesh::err << "WARNING: InfFE<Dim,T_radial,T_map>::shape(...) does _not_" << std::endl
215  << " return the correct trial function! Use " << std::endl
216  << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!"
217  << std::endl;
218  _warned_for_shape = true;
219  }
220 #endif
221 
222  const Order o_radial (fet.radial_order);
223  const Real v (p(Dim-1));
224  std::unique_ptr<const Elem> base_el (inf_elem->build_side_ptr(0));
225 
226  unsigned int i_base, i_radial;
227  compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial);
228 
229  if (Dim > 1)
230  return FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p)
231  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
232  * InfFERadial::decay(Dim,v);
233  else
234  return InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
235  * InfFERadial::decay(Dim,v);
236 }
237 
238 
239 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
241  const ElemType inf_elem_type,
242  const unsigned int i,
243  const unsigned int j,
244  const Point & p)
245 {
246  libmesh_assert_not_equal_to (Dim, 0);
247  libmesh_assert_greater (Dim,j);
248 #ifdef DEBUG
249  // this makes only sense when used for mapping
250  if ((T_radial != INFINITE_MAP) && !_warned_for_dshape)
251  {
252  libMesh::err << "WARNING: InfFE<Dim,T_radial,T_map>::shape_deriv(...) does _not_" << std::endl
253  << " return the correct trial function gradients! Use " << std::endl
254  << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!"
255  << std::endl;
256  _warned_for_dshape = true;
257  }
258 #endif
259 
260  const ElemType base_et (InfFEBase::get_elem_type(inf_elem_type));
261  const Order o_radial (fe_t.radial_order);
262  const Real v (p(Dim-1));
263 
264  unsigned int i_base, i_radial;
265  compute_shape_indices(fe_t, inf_elem_type, i, i_base, i_radial);
266 
267  if (j== Dim -1)
268  {
269  Real RadialDeriv = InfFE<Dim,T_radial,T_map>::eval_deriv(v, o_radial, i_radial)
270  * InfFERadial::decay(Dim,v)
271  + InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
273 
274  return FEInterface::shape(Dim-1, fe_t, base_et, i_base, p)*RadialDeriv;
275  }
276 
277  return FEInterface::shape_deriv(Dim-1, fe_t, base_et, i_base, j, p)
278  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
279  * InfFERadial::decay(Dim,v);
280 }
281 
282 
283 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
285  const Elem * inf_elem,
286  const unsigned int i,
287  const unsigned int j,
288  const Point & p)
289 {
290  libmesh_assert_not_equal_to (Dim, 0);
291  libmesh_assert_greater (Dim,j);
292 #ifdef DEBUG
293  // this makes only sense when used for mapping
294  if ((T_radial != INFINITE_MAP) && !_warned_for_dshape)
295  {
296  libMesh::err << "WARNING: InfFE<Dim,T_radial,T_map>::shape_deriv(...) does _not_" << std::endl
297  << " return the correct trial function gradients! Use " << std::endl
298  << " InfFE<Dim,T_radial,T_map>::compute_data(..) instead!"
299  << std::endl;
300  _warned_for_dshape = true;
301  }
302 #endif
303  const Order o_radial (fe_t.radial_order);
304  const Real v (p(Dim-1));
305 
306  std::unique_ptr<const Elem> base_el (inf_elem->build_side_ptr(0));
307 
308  unsigned int i_base, i_radial;
309 
310  if ((-1. > v ) || (v > 1.))
311  {
312  //TODO: This is for debugging. We should never come here.
313  // Therefore we can do very useless things then:
314  i_base=0;
315  }
316  compute_shape_indices(fe_t, inf_elem->type(), i, i_base, i_radial);
317 
318  if (j== Dim -1)
319  {
320  Real RadialDeriv = InfFE<Dim,T_radial,T_map>::eval_deriv(v, o_radial, i_radial)
321  * InfFERadial::decay(Dim,v)
322  + InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
324 
325  return FEInterface::shape(Dim-1, fe_t, base_el.get(), i_base, p)*RadialDeriv;
326  }
327  return FEInterface::shape_deriv(Dim-1, fe_t, base_el.get(), i_base, j, p)
328  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
329  * InfFERadial::decay(Dim,v);
330 }
331 
332 
333 
334 
335 
336 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
338  const Elem * inf_elem,
340 {
341  libmesh_assert(inf_elem);
342  libmesh_assert_not_equal_to (Dim, 0);
343 
344  const Order o_radial (fet.radial_order);
345  const Order radial_mapping_order (InfFERadial::mapping_order());
346  const Point & p (data.p);
347  const Real v (p(Dim-1));
348  std::unique_ptr<const Elem> base_el (inf_elem->build_side_ptr(0));
349 
350  /*
351  * compute \p interpolated_dist containing the mapping-interpolated
352  * distance of the base point to the origin. This is the same
353  * for all shape functions. Set \p interpolated_dist to 0, it
354  * is added to.
355  */
356  Real interpolated_dist = 0.;
357  switch (Dim)
358  {
359  case 1:
360  {
361  libmesh_assert_equal_to (inf_elem->type(), INFEDGE2);
362  interpolated_dist = Point(inf_elem->point(0) - inf_elem->point(1)).norm();
363  break;
364  }
365 
366  case 2:
367  {
368  const unsigned int n_base_nodes = base_el->n_nodes();
369 
370  const Point origin = inf_elem->origin();
371  const Order base_mapping_order (base_el->default_order());
372  const ElemType base_mapping_elem_type (base_el->type());
373 
374  // interpolate the base nodes' distances
375  for (unsigned int n=0; n<n_base_nodes; n++)
376  interpolated_dist += Point(base_el->point(n) - origin).norm()
377  * FE<1,LAGRANGE>::shape (base_mapping_elem_type, base_mapping_order, n, p);
378  break;
379  }
380 
381  case 3:
382  {
383  const unsigned int n_base_nodes = base_el->n_nodes();
384 
385  const Point origin = inf_elem->origin();
386  const Order base_mapping_order (base_el->default_order());
387  const ElemType base_mapping_elem_type (base_el->type());
388 
389  // interpolate the base nodes' distances
390  for (unsigned int n=0; n<n_base_nodes; n++)
391  interpolated_dist += Point(base_el->point(n) - origin).norm()
392  * FE<2,LAGRANGE>::shape (base_mapping_elem_type, base_mapping_order, n, p);
393  break;
394  }
395 
396  default:
397  libmesh_error_msg("Unknown Dim = " << Dim);
398  }
399 
400 
401  const Real speed = data.speed;
402 
403  //TODO: I find it inconvenient to have a quantity phase which is phase/speed.
404  // But it might be better than redefining a quantities meaning.
405  data.phase = interpolated_dist /* together with next line: */
406  * InfFE<Dim,INFINITE_MAP,T_map>::eval(v, radial_mapping_order, 1)/speed; /* phase(s,t,v)/c */
407 
408  // We assume time-harmonic behavior in this function!
409 
410 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
411  // the wave number
412  const Number wavenumber = 2. * libMesh::pi * data.frequency / speed;
413 
414  // the exponent for time-harmonic behavior
415  // \note: this form is much less general than the implementation of dphase, which can be easily extended to
416  // other forms than e^{i kr}.
417  const Number exponent = imaginary /* imaginary unit */
418  * wavenumber /* k (can be complex) */
419  * data.phase*speed;
420 
421  const Number time_harmonic = exp(exponent); /* e^(i*k*phase(s,t,v)) */
422 #else
423  const Number time_harmonic = 1;
424 #endif //LIBMESH_USE_COMPLEX_NUMBERS
425 
426  /*
427  * compute \p shape for all dof in the element
428  */
429  if (Dim > 1)
430  {
431  const unsigned int n_dof = n_dofs (fet, inf_elem->type());
432  data.shape.resize(n_dof);
433  if (data.need_derivative())
434  {
435  data.dshape.resize(n_dof);
436  data.local_transform.resize(Dim);
437 
438  for (unsigned int d=0; d<Dim; d++)
439  data.local_transform[d].resize(Dim);
440 
441  // compute the reference->physical map at the point \p p.
442  // Use another fe_map to avoid interference with \p this->_fe_map
443  // which is initialized at the quadrature points...
444  UniquePtr<FEBase> fe (FEBase::build_InfFE(Dim, fet));
445  std::vector<Point> pt(1);
446  pt[0]=p;
447  fe->get_dphideta(); // to compute the map
448  fe->reinit(inf_elem, &pt);
449 
450  // compute the reference->physical map.
451  data.local_transform[0][0] = fe->get_dxidx()[0];
452  data.local_transform[1][0] = fe->get_detadx()[0];
453  data.local_transform[1][1] = fe->get_detady()[0];
454  data.local_transform[0][1] = fe->get_dxidy()[0];
455  if (Dim > 2)
456  {
457  data.local_transform[2][0] = fe->get_dzetadx()[0];
458  data.local_transform[2][1] = fe->get_dzetady()[0];
459  data.local_transform[2][2] = fe->get_dzetadz()[0];
460  data.local_transform[1][2] = fe->get_detadz()[0];
461  data.local_transform[0][2] = fe->get_dxidz()[0];
462  }
463  } // endif data.need_derivative()
464 
465  for (unsigned int i=0; i<n_dof; i++)
466  {
467  // compute base and radial shape indices
468  unsigned int i_base, i_radial;
469  compute_shape_indices(fet, inf_elem->type(), i, i_base, i_radial);
470 
471  data.shape[i] = (InfFERadial::decay(Dim,v) /* (1.-v)/2. in 3D */
472  * FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) /* S_n(s,t) */
473  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)) /* L_n(v) */
474  * time_harmonic; /* e^(i*k*phase(s,t,v) */
475 
476  // use differentiation of the above equation
477  if (data.need_derivative())
478  {
479  data.dshape[i](0) = (InfFERadial::decay(Dim,v)
480  * FEInterface::shape_deriv(Dim-1, fet, base_el.get(), i_base, 0, p)
481  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial))
482  * time_harmonic;
483 
484  if (Dim > 2)
485  {
486  data.dshape[i](1) = (InfFERadial::decay(Dim,v)
487  * FEInterface::shape_deriv(Dim-1, fet, base_el.get(), i_base, 1, p)
488  * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial))
489  * time_harmonic;
490 
491  }
492  data.dshape[i](Dim-1) = (InfFERadial::decay_deriv(v) * InfFE<Dim,T_radial,T_map>::eval(v, o_radial, i_radial)
493  +InfFERadial::decay(Dim,v) * InfFE<Dim,T_radial,T_map>::eval_deriv(v, o_radial, i_radial))
494  * FEInterface::shape(Dim-1, fet, base_el.get(), i_base, p) * time_harmonic;
495 
496 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
497  // derivative of time_harmonic (works for harmonic behavior only):
498  data.dshape[i](Dim-1)+= data.shape[i]*imaginary*wavenumber
499  *interpolated_dist*InfFE<Dim,INFINITE_MAP,T_map>::eval_deriv(v, radial_mapping_order, 1);
500 
501 #else
502  /*
503  * The gradient in infinite elements is dominated by the contribution due to the oscillating phase.
504  * Since this term is imaginary, I think there is no means to look at it without having complex numbers.
505  */
506  libmesh_not_implemented();
507  // Maybe we can solve it with a warning as well, but I think one really should not do this...
508 #endif
509  }
510  }
511  }
512 
513  else
514  libmesh_error_msg("compute_data() for 1-dimensional InfFE not implemented.");
515 }
516 
517 
518 
519 
520 
521 
522 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
524  const unsigned int outer_node_index,
525  unsigned int & base_node,
526  unsigned int & radial_node)
527 {
528  switch (inf_elem_type)
529  {
530  case INFEDGE2:
531  {
532  libmesh_assert_less (outer_node_index, 2);
533  base_node = 0;
534  radial_node = outer_node_index;
535  return;
536  }
537 
538 
539  // linear base approximation, easy to determine
540  case INFQUAD4:
541  {
542  libmesh_assert_less (outer_node_index, 4);
543  base_node = outer_node_index % 2;
544  radial_node = outer_node_index / 2;
545  return;
546  }
547 
548  case INFPRISM6:
549  {
550  libmesh_assert_less (outer_node_index, 6);
551  base_node = outer_node_index % 3;
552  radial_node = outer_node_index / 3;
553  return;
554  }
555 
556  case INFHEX8:
557  {
558  libmesh_assert_less (outer_node_index, 8);
559  base_node = outer_node_index % 4;
560  radial_node = outer_node_index / 4;
561  return;
562  }
563 
564 
565  // higher order base approximation, more work necessary
566  case INFQUAD6:
567  {
568  switch (outer_node_index)
569  {
570  case 0:
571  case 1:
572  {
573  radial_node = 0;
574  base_node = outer_node_index;
575  return;
576  }
577 
578  case 2:
579  case 3:
580  {
581  radial_node = 1;
582  base_node = outer_node_index-2;
583  return;
584  }
585 
586  case 4:
587  {
588  radial_node = 0;
589  base_node = 2;
590  return;
591  }
592 
593  case 5:
594  {
595  radial_node = 1;
596  base_node = 2;
597  return;
598  }
599 
600  default:
601  libmesh_error_msg("Unrecognized outer_node_index = " << outer_node_index);
602  }
603  }
604 
605 
606  case INFHEX16:
607  case INFHEX18:
608  {
609  switch (outer_node_index)
610  {
611  case 0:
612  case 1:
613  case 2:
614  case 3:
615  {
616  radial_node = 0;
617  base_node = outer_node_index;
618  return;
619  }
620 
621  case 4:
622  case 5:
623  case 6:
624  case 7:
625  {
626  radial_node = 1;
627  base_node = outer_node_index-4;
628  return;
629  }
630 
631  case 8:
632  case 9:
633  case 10:
634  case 11:
635  {
636  radial_node = 0;
637  base_node = outer_node_index-4;
638  return;
639  }
640 
641  case 12:
642  case 13:
643  case 14:
644  case 15:
645  {
646  radial_node = 1;
647  base_node = outer_node_index-8;
648  return;
649  }
650 
651  case 16:
652  {
653  libmesh_assert_equal_to (inf_elem_type, INFHEX18);
654  radial_node = 0;
655  base_node = 8;
656  return;
657  }
658 
659  case 17:
660  {
661  libmesh_assert_equal_to (inf_elem_type, INFHEX18);
662  radial_node = 1;
663  base_node = 8;
664  return;
665  }
666 
667  default:
668  libmesh_error_msg("Unrecognized outer_node_index = " << outer_node_index);
669  }
670  }
671 
672 
673  case INFPRISM12:
674  {
675  switch (outer_node_index)
676  {
677  case 0:
678  case 1:
679  case 2:
680  {
681  radial_node = 0;
682  base_node = outer_node_index;
683  return;
684  }
685 
686  case 3:
687  case 4:
688  case 5:
689  {
690  radial_node = 1;
691  base_node = outer_node_index-3;
692  return;
693  }
694 
695  case 6:
696  case 7:
697  case 8:
698  {
699  radial_node = 0;
700  base_node = outer_node_index-3;
701  return;
702  }
703 
704  case 9:
705  case 10:
706  case 11:
707  {
708  radial_node = 1;
709  base_node = outer_node_index-6;
710  return;
711  }
712 
713  default:
714  libmesh_error_msg("Unrecognized outer_node_index = " << outer_node_index);
715  }
716  }
717 
718 
719  default:
720  libmesh_error_msg("ERROR: Bad infinite element type=" << inf_elem_type << ", node=" << outer_node_index);
721  }
722 }
723 
724 
725 
726 
727 
728 
729 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
731  const unsigned int outer_node_index,
732  unsigned int & base_node,
733  unsigned int & radial_node)
734 {
735  libmesh_assert_not_equal_to (inf_elem_type, INVALID_ELEM);
736 
737  static std::vector<unsigned int> _static_base_node_index;
738  static std::vector<unsigned int> _static_radial_node_index;
739 
740  /*
741  * fast counterpart to compute_node_indices(), uses local static buffers
742  * to store the index maps. The class member
743  * \p _compute_node_indices_fast_current_elem_type remembers
744  * the current element type.
745  *
746  * Note that there exist non-static members storing the
747  * same data. However, you never know what element type
748  * is currently used by the \p InfFE object, and what
749  * request is currently directed to the static \p InfFE
750  * members (which use \p compute_node_indices_fast()).
751  * So separate these.
752  *
753  * check whether the work for this elemtype has already
754  * been done. If so, use this index. Otherwise, refresh
755  * the buffer to this element type.
756  */
757  if (inf_elem_type==_compute_node_indices_fast_current_elem_type)
758  {
759  base_node = _static_base_node_index [outer_node_index];
760  radial_node = _static_radial_node_index[outer_node_index];
761  return;
762  }
763  else
764  {
765  // store the map for _all_ nodes for this element type
766  _compute_node_indices_fast_current_elem_type = inf_elem_type;
767 
768  unsigned int n_nodes = libMesh::invalid_uint;
769 
770  switch (inf_elem_type)
771  {
772  case INFEDGE2:
773  {
774  n_nodes = 2;
775  break;
776  }
777  case INFQUAD4:
778  {
779  n_nodes = 4;
780  break;
781  }
782  case INFQUAD6:
783  {
784  n_nodes = 6;
785  break;
786  }
787  case INFHEX8:
788  {
789  n_nodes = 8;
790  break;
791  }
792  case INFHEX16:
793  {
794  n_nodes = 16;
795  break;
796  }
797  case INFHEX18:
798  {
799  n_nodes = 18;
800  break;
801  }
802  case INFPRISM6:
803  {
804  n_nodes = 6;
805  break;
806  }
807  case INFPRISM12:
808  {
809  n_nodes = 12;
810  break;
811  }
812  default:
813  libmesh_error_msg("ERROR: Bad infinite element type=" << inf_elem_type << ", node=" << outer_node_index);
814  }
815 
816 
817  _static_base_node_index.resize (n_nodes);
818  _static_radial_node_index.resize(n_nodes);
819 
820  for (unsigned int n=0; n<n_nodes; n++)
821  compute_node_indices (inf_elem_type,
822  n,
823  _static_base_node_index [outer_node_index],
824  _static_radial_node_index[outer_node_index]);
825 
826  // and return for the specified node
827  base_node = _static_base_node_index [outer_node_index];
828  radial_node = _static_radial_node_index[outer_node_index];
829  return;
830  }
831 }
832 
833 
834 
835 
836 
837 
838 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
840  const ElemType inf_elem_type,
841  const unsigned int i,
842  unsigned int & base_shape,
843  unsigned int & radial_shape)
844 {
845 
846  /*
847  * An example is provided: the numbers in comments refer to
848  * a fictitious InfHex18. The numbers are chosen as exemplary
849  * values. There is currently no base approximation that
850  * requires this many dof's at nodes, sides, faces and in the element.
851  *
852  * the order of the shape functions is heavily related with the
853  * order the dofs are assigned in \p DofMap::distributed_dofs().
854  * Due to the infinite elements with higher-order base approximation,
855  * some more effort is necessary.
856  *
857  * numbering scheme:
858  * 1. all vertices in the base, assign node->n_comp() dofs to each vertex
859  * 2. all vertices further out: innermost loop: radial shapes,
860  * then the base approximation shapes
861  * 3. all side nodes in the base, assign node->n_comp() dofs to each side node
862  * 4. all side nodes further out: innermost loop: radial shapes,
863  * then the base approximation shapes
864  * 5. (all) face nodes in the base, assign node->n_comp() dofs to each face node
865  * 6. (all) face nodes further out: innermost loop: radial shapes,
866  * then the base approximation shapes
867  * 7. element-associated dof in the base
868  * 8. element-associated dof further out
869  */
870 
871  const unsigned int radial_order = static_cast<unsigned int>(fet.radial_order.get_order()); // 4
872  const unsigned int radial_order_p_one = radial_order+1; // 5
873 
874  const ElemType base_elem_type (InfFEBase::get_elem_type(inf_elem_type)); // QUAD9
875 
876  // assume that the number of dof is the same for all vertices
877  unsigned int n_base_vertices = libMesh::invalid_uint; // 4
878  const unsigned int n_base_vertex_dof = FEInterface::n_dofs_at_node (Dim-1, fet, base_elem_type, 0);// 2
879 
880  unsigned int n_base_side_nodes = libMesh::invalid_uint; // 4
881  unsigned int n_base_side_dof = libMesh::invalid_uint; // 3
882 
883  unsigned int n_base_face_nodes = libMesh::invalid_uint; // 1
884  unsigned int n_base_face_dof = libMesh::invalid_uint; // 5
885 
886  const unsigned int n_base_elem_dof = FEInterface::n_dofs_per_elem (Dim-1, fet, base_elem_type);// 9
887 
888 
889  switch (inf_elem_type)
890  {
891  case INFEDGE2:
892  {
893  n_base_vertices = 1;
894  n_base_side_nodes = 0;
895  n_base_face_nodes = 0;
896  n_base_side_dof = 0;
897  n_base_face_dof = 0;
898  break;
899  }
900 
901  case INFQUAD4:
902  {
903  n_base_vertices = 2;
904  n_base_side_nodes = 0;
905  n_base_face_nodes = 0;
906  n_base_side_dof = 0;
907  n_base_face_dof = 0;
908  break;
909  }
910 
911  case INFQUAD6:
912  {
913  n_base_vertices = 2;
914  n_base_side_nodes = 1;
915  n_base_face_nodes = 0;
916  n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices);
917  n_base_face_dof = 0;
918  break;
919  }
920 
921  case INFHEX8:
922  {
923  n_base_vertices = 4;
924  n_base_side_nodes = 0;
925  n_base_face_nodes = 0;
926  n_base_side_dof = 0;
927  n_base_face_dof = 0;
928  break;
929  }
930 
931  case INFHEX16:
932  {
933  n_base_vertices = 4;
934  n_base_side_nodes = 4;
935  n_base_face_nodes = 0;
936  n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices);
937  n_base_face_dof = 0;
938  break;
939  }
940 
941  case INFHEX18:
942  {
943  n_base_vertices = 4;
944  n_base_side_nodes = 4;
945  n_base_face_nodes = 1;
946  n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices);
947  n_base_face_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, 8);
948  break;
949  }
950 
951 
952  case INFPRISM6:
953  {
954  n_base_vertices = 3;
955  n_base_side_nodes = 0;
956  n_base_face_nodes = 0;
957  n_base_side_dof = 0;
958  n_base_face_dof = 0;
959  break;
960  }
961 
962  case INFPRISM12:
963  {
964  n_base_vertices = 3;
965  n_base_side_nodes = 3;
966  n_base_face_nodes = 0;
967  n_base_side_dof = FEInterface::n_dofs_at_node (Dim-1, fet,base_elem_type, n_base_vertices);
968  n_base_face_dof = 0;
969  break;
970  }
971 
972  default:
973  libmesh_error_msg("Unrecognized inf_elem_type = " << inf_elem_type);
974  }
975 
976 
977  {
978  // these are the limits describing the intervals where the shape function lies
979  const unsigned int n_dof_at_base_vertices = n_base_vertices*n_base_vertex_dof; // 8
980  const unsigned int n_dof_at_all_vertices = n_dof_at_base_vertices*radial_order_p_one; // 40
981 
982  const unsigned int n_dof_at_base_sides = n_base_side_nodes*n_base_side_dof; // 12
983  const unsigned int n_dof_at_all_sides = n_dof_at_base_sides*radial_order_p_one; // 60
984 
985  const unsigned int n_dof_at_base_face = n_base_face_nodes*n_base_face_dof; // 5
986  const unsigned int n_dof_at_all_faces = n_dof_at_base_face*radial_order_p_one; // 25
987 
988 
989  // start locating the shape function
990  if (i < n_dof_at_base_vertices) // range of i: 0..7
991  {
992  // belongs to vertex in the base
993  radial_shape = 0;
994  base_shape = i;
995  }
996 
997  else if (i < n_dof_at_all_vertices) // range of i: 8..39
998  {
999  /* belongs to vertex in the outer shell
1000  *
1001  * subtract the number of dof already counted,
1002  * so that i_offset contains only the offset for the base
1003  */
1004  const unsigned int i_offset = i - n_dof_at_base_vertices; // 0..31
1005 
1006  // first the radial dof are counted, then the base dof
1007  radial_shape = (i_offset % radial_order) + 1;
1008  base_shape = i_offset / radial_order;
1009  }
1010 
1011  else if (i < n_dof_at_all_vertices+n_dof_at_base_sides) // range of i: 40..51
1012  {
1013  // belongs to base, is a side node
1014  radial_shape = 0;
1015  base_shape = i - radial_order * n_dof_at_base_vertices; // 8..19
1016  }
1017 
1018  else if (i < n_dof_at_all_vertices+n_dof_at_all_sides) // range of i: 52..99
1019  {
1020  // belongs to side node in the outer shell
1021  const unsigned int i_offset = i - (n_dof_at_all_vertices
1022  + n_dof_at_base_sides); // 0..47
1023  radial_shape = (i_offset % radial_order) + 1;
1024  base_shape = (i_offset / radial_order) + n_dof_at_base_vertices;
1025  }
1026 
1027  else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_base_face) // range of i: 100..104
1028  {
1029  // belongs to the node in the base face
1030  radial_shape = 0;
1031  base_shape = i - radial_order*(n_dof_at_base_vertices
1032  + n_dof_at_base_sides); // 20..24
1033  }
1034 
1035  else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_all_faces) // range of i: 105..124
1036  {
1037  // belongs to the node in the outer face
1038  const unsigned int i_offset = i - (n_dof_at_all_vertices
1039  + n_dof_at_all_sides
1040  + n_dof_at_base_face); // 0..19
1041  radial_shape = (i_offset % radial_order) + 1;
1042  base_shape = (i_offset / radial_order) + n_dof_at_base_vertices + n_dof_at_base_sides;
1043  }
1044 
1045  else if (i < n_dof_at_all_vertices+n_dof_at_all_sides+n_dof_at_all_faces+n_base_elem_dof) // range of i: 125..133
1046  {
1047  // belongs to the base and is an element associated shape
1048  radial_shape = 0;
1049  base_shape = i - (n_dof_at_all_vertices
1050  + n_dof_at_all_sides
1051  + n_dof_at_all_faces); // 0..8
1052  }
1053 
1054  else // range of i: 134..169
1055  {
1056  libmesh_assert_less (i, n_dofs(fet, inf_elem_type));
1057  // belongs to the outer shell and is an element associated shape
1058  const unsigned int i_offset = i - (n_dof_at_all_vertices
1059  + n_dof_at_all_sides
1060  + n_dof_at_all_faces
1061  + n_base_elem_dof); // 0..19
1062  radial_shape = (i_offset % radial_order) + 1;
1063  base_shape = (i_offset / radial_order) + n_dof_at_base_vertices + n_dof_at_base_sides + n_dof_at_base_face;
1064  }
1065  }
1066 
1067  return;
1068 }
1069 
1070 
1071 
1072 //--------------------------------------------------------------
1073 // Explicit instantiations
1074 //#include "libmesh/inf_fe_instantiate_1D.h"
1075 //#include "libmesh/inf_fe_instantiate_2D.h"
1076 //#include "libmesh/inf_fe_instantiate_3D.h"
1077 
1078 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, unsigned int, n_dofs(const FEType &, const ElemType));
1079 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, unsigned int, n_dofs(const FEType &, const ElemType));
1080 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, unsigned int, n_dofs(const FEType &, const ElemType));
1081 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, unsigned int, n_dofs_per_elem(const FEType &, const ElemType));
1082 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, unsigned int, n_dofs_per_elem(const FEType &, const ElemType));
1083 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, unsigned int, n_dofs_per_elem(const FEType &, const ElemType));
1084 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, unsigned int, n_dofs_at_node(const FEType &, const ElemType, const unsigned int));
1085 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, unsigned int, n_dofs_at_node(const FEType &, const ElemType, const unsigned int));
1086 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, unsigned int, n_dofs_at_node(const FEType &, const ElemType, const unsigned int));
1087 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, void, compute_shape_indices(const FEType &, const ElemType, const unsigned int, unsigned int &, unsigned int &));
1088 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, void, compute_shape_indices(const FEType &, const ElemType, const unsigned int, unsigned int &, unsigned int &));
1089 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, void, compute_shape_indices(const FEType &, const ElemType, const unsigned int, unsigned int &, unsigned int &));
1090 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, void, compute_node_indices(const ElemType, const unsigned int, unsigned int &, unsigned int &));
1091 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, void, compute_node_indices(const ElemType, const unsigned int, unsigned int &, unsigned int &));
1092 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, void, compute_node_indices(const ElemType, const unsigned int, unsigned int &, unsigned int &));
1093 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, Real, shape(const FEType &, const Elem *, const unsigned int, const Point & p));
1094 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, Real, shape(const FEType &, const Elem *, const unsigned int, const Point & p));
1095 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, Real, shape(const FEType &, const Elem *, const unsigned int, const Point & p));
1096 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, Real, shape(const FEType &, const ElemType, const unsigned int, const Point &));
1097 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, Real, shape(const FEType &, const ElemType, const unsigned int, const Point &));
1098 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, Real, shape(const FEType &, const ElemType, const unsigned int, const Point &));
1099 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, Real, shape_deriv(const FEType &, const Elem *, const unsigned int, const unsigned int, const Point &));
1100 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, Real, shape_deriv(const FEType &, const Elem *, const unsigned int, const unsigned int, const Point &));
1101 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, Real, shape_deriv(const FEType &, const Elem *, const unsigned int, const unsigned int, const Point &));
1102 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, Real, shape_deriv(const FEType &, const ElemType, const unsigned int, const unsigned int, const Point &));
1103 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, Real, shape_deriv(const FEType &, const ElemType, const unsigned int, const unsigned int, const Point &));
1104 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, Real, shape_deriv(const FEType &, const ElemType, const unsigned int, const unsigned int, const Point &));
1105 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, void, compute_data(const FEType &, const Elem *, FEComputeData &));
1106 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, void, compute_data(const FEType &, const Elem *, FEComputeData &));
1107 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, void, compute_data(const FEType &, const Elem *, FEComputeData &));
1108 INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, void, nodal_soln(const FEType &, const Elem *, const std::vector<Number> &, std::vector<Number> &));
1109 INSTANTIATE_INF_FE_MBRF(2, CARTESIAN, void, nodal_soln(const FEType &, const Elem *, const std::vector<Number> &, std::vector<Number> &));
1110 INSTANTIATE_INF_FE_MBRF(3, CARTESIAN, void, nodal_soln(const FEType &, const Elem *, const std::vector<Number> &, std::vector<Number> &));
1111 
1112 } // namespace libMesh
1113 
1114 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
libMesh::InfFERadial::decay
static Real decay(const unsigned int dim, const Real v)
Definition: inf_fe.h:850
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::pi
const Real pi
.
Definition: libmesh.h:237
libMesh::OrderWrapper::get_order
int get_order() const
Explicitly request the order as an int.
Definition: fe_type.h:77
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::InfFERadial::mapping_order
static Order mapping_order()
Definition: inf_fe.h:93
libMesh::FEComputeData
class FEComputeData hides arbitrary data to be passed to and from children of FEBase through the FEIn...
Definition: fe_compute_data.h:51
libMesh::Elem::build_side_ptr
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i, bool proxy=true)=0
libMesh::INFINITE_MAP
Definition: enum_fe_family.h:48
libMesh::INFHEX8
Definition: enum_elem_type.h:60
libMesh::InfFE
A specific instantiation of the FEBase class.
Definition: fe.h:40
libMesh::INFQUAD4
Definition: enum_elem_type.h:58
libMesh::InfFE::eval_deriv
static Real eval_deriv(Real v, Order o_radial, unsigned int i)
libMesh::FEInterface::n_dofs_at_node
static unsigned int n_dofs_at_node(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int n)
Definition: fe_interface.C:503
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::InfFE::n_dofs
static unsigned int n_dofs(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:58
libMesh::InfFE::nodal_soln
static void nodal_soln(const FEType &fet, const Elem *elem, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Usually, this method would build the nodal soln from the element soln.
Definition: inf_fe_static.C:123
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::INFEDGE2
Definition: enum_elem_type.h:57
libMesh::FEGenericBase::build_InfFE
static std::unique_ptr< FEGenericBase > build_InfFE(const unsigned int dim, const FEType &type)
Builds a specific infinite element type.
libMesh::Elem::point
const Point & point(const unsigned int i) const
Definition: elem.h:1955
libMesh::CARTESIAN
Definition: enum_inf_map_type.h:35
libMesh::INSTANTIATE_INF_FE_MBRF
INSTANTIATE_INF_FE_MBRF(1, CARTESIAN, void, reinit(const Elem *, const unsigned int, const Real, const std::vector< Point > *const, const std::vector< Real > *const))
libMesh::InfFERadial::decay_deriv
static Real decay_deriv(const Real)
Definition: inf_fe.h:76
libMesh::INFPRISM6
Definition: enum_elem_type.h:63
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::imaginary
const Number imaginary
The imaginary unit, .
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::UniquePtr
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:33
libMesh::FEInterface::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:534
libMesh::InfFE::compute_node_indices
static void compute_node_indices(const ElemType inf_elem_type, const unsigned int outer_node_index, unsigned int &base_node, unsigned int &radial_node)
Computes the indices in the base base_node and in radial direction radial_node (either 0 or 1) associ...
Definition: inf_fe_static.C:523
libMesh::InfFE::_warned_for_dshape
static bool _warned_for_dshape
Definition: inf_fe.h:831
libMesh::INFHEX18
Definition: enum_elem_type.h:62
libMesh::InfFERadial::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const Order o_radial)
Definition: inf_fe.h:132
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
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
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::INFHEX16
Definition: enum_elem_type.h:61
libMesh::InfFE::_warned_for_nodal_soln
static bool _warned_for_nodal_soln
static members that are used to issue warning messages only once.
Definition: inf_fe.h:829
libMesh::InfFERadial::n_dofs
static unsigned int n_dofs(const Order o_radial)
Definition: inf_fe.h:109
libMesh::INFPRISM12
Definition: enum_elem_type.h:64
libMesh::INVALID_ELEM
Definition: enum_elem_type.h:75
libMesh::InfFE::shape_deriv
static Real shape_deriv(const FEType &fet, const Elem *inf_elem, const unsigned int i, const unsigned int j, const Point &p)
Definition: inf_fe_static.C:284
libMesh::FEInterface::shape_deriv
static Real shape_deriv(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const unsigned int j, const Point &p)
Definition: fe_interface.C:845
libMesh::FEType
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:178
libMesh::INFQUAD6
Definition: enum_elem_type.h:59
libMesh::InfFE::shape
static Real shape(const FEType &fet, const ElemType t, const unsigned int i, const Point &p)
Definition: inf_fe_static.C:158
libMesh::InfFERadial::n_dofs_at_node
static unsigned int n_dofs_at_node(const Order o_radial, const unsigned int n_onion)
Definition: inf_fe_base_radial.C:119
libMesh::InfFE::_compute_node_indices_fast_current_elem_type
static ElemType _compute_node_indices_fast_current_elem_type
When compute_node_indices_fast() is used, this static variable remembers the element type for which t...
Definition: inf_fe.h:821
libMesh::InfFE::compute_shape_indices
static void compute_shape_indices(const FEType &fet, const ElemType inf_elem_type, const unsigned int i, unsigned int &base_shape, unsigned int &radial_shape)
Computes the indices of shape functions in the base base_shape and in radial direction radial_shape (...
Definition: inf_fe_static.C:839
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
std::norm
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
data
IterBase * data
Ideally this private member data should have protected access.
Definition: variant_filter_iterator.h:337
libMesh::Elem::origin
virtual Point origin() const
Definition: elem.h:1593
libMesh::err
OStreamProxy err
libMesh::InfFE::eval
static Real eval(Real v, Order o_radial, unsigned int i)
libMesh::InfFE::n_dofs_at_node
static unsigned int n_dofs_at_node(const FEType &fet, const ElemType inf_elem_type, const unsigned int n)
Definition: inf_fe_static.C:76
libMesh::InfFE::compute_data
static void compute_data(const FEType &fe_t, const Elem *inf_elem, FEComputeData &data)
Generalized version of shape(), takes an Elem *.
Definition: inf_fe_static.C:337
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::InfFE::compute_node_indices_fast
static void compute_node_indices_fast(const ElemType inf_elem_type, const unsigned int outer_node_index, unsigned int &base_node, unsigned int &radial_node)
Does the same as compute_node_indices(), but stores the maps for the current element type.
Definition: inf_fe_static.C:730
libMesh::InfFE::_warned_for_shape
static bool _warned_for_shape
Definition: inf_fe.h:830
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::FEType::radial_order
OrderWrapper radial_order
The approximation order in radial direction of the infinite element.
Definition: fe_type.h:236
libMesh::InfFEBase::get_elem_type
static ElemType get_elem_type(const ElemType type)
Definition: inf_fe_base_radial.C:49
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
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
libMesh::InfFE::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const FEType &fet, const ElemType inf_elem_type)
Definition: inf_fe_static.C:105