libMesh
inf_fe.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2024 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/libmesh_config.h"
22 
23 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
24 
25 #include "libmesh/inf_fe.h"
26 #include "libmesh/fe.h"
27 #include "libmesh/quadrature_gauss.h"
28 #include "libmesh/elem.h"
29 #include "libmesh/libmesh_logging.h"
30 #include "libmesh/int_range.h"
31 #include "libmesh/type_tensor.h"
32 #include "libmesh/fe_interface.h"
33 
34 #include <memory>
35 
36 namespace libMesh
37 {
38 
39 
40 
41 // Constructor
42 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
44  FEBase (Dim, fet),
45 
46  calculate_map_scaled(false),
47  calculate_phi_scaled(false),
48  calculate_dphi_scaled(false),
49  calculate_xyz(false),
50  calculate_jxw(false),
51  _n_total_approx_sf (0),
52  _n_total_qp (0),
53 
54  // initialize the current_fe_type to all the same
55  // values as \p fet (since the FE families and coordinate
56  // map type should not change), but use an invalid order
57  // for the radial part (since this is the only order
58  // that may change!).
59  // the data structures like \p phi etc are not initialized
60  // through the constructor, but through reinit()
61  current_fe_type (FEType(fet.order,
62  fet.family,
64  fet.radial_family,
65  fet.inf_map))
66 
67 {
68  // Sanity checks
69  libmesh_assert_equal_to (T_radial, fe_type.radial_family);
70  libmesh_assert_equal_to (T_map, fe_type.inf_map);
71 
72  // build the base_fe object
73  if (Dim != 1)
74  base_fe = FEBase::build(Dim-1, fet);
75 }
76 
77 
78 
79 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
81 {
82  libmesh_assert(q);
83  libmesh_assert(base_fe);
84 
85  const Order base_int_order = q->get_order();
86  const Order radial_int_order = static_cast<Order>(2 * (static_cast<unsigned int>(fe_type.radial_order.get_order()) + 1) +2);
87  const unsigned int qrule_dim = q->get_dim();
88 
89  if (Dim != 1)
90  {
91  // build a Dim-1 quadrature rule of the type that we received
92  base_qrule = QBase::build(q->type(), qrule_dim-1, base_int_order);
93  base_fe->attach_quadrature_rule(base_qrule.get());
94  }
95 
96  // in radial direction, always use Gauss quadrature
97  radial_qrule = std::make_unique<QGauss>(1, radial_int_order);
98 
99  // Maybe helpful to store the QBase *
100  // with which we initialized our own quadrature rules.
101  // Used e.g. in \p InfFE::reinit(elem,side)
102  qrule = q;
103 }
104 
105 
106 
107 
108 
109 template <unsigned int Dim, FEFamily T_radial, InfMapType T_base>
111 {
112  base_elem = InfFEBase::build_elem(inf_elem);
113 }
114 
115 
116 
117 
118 
119 
120 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
122  const std::vector<Point> * const pts,
123  const std::vector<Real> * const weights)
124 {
125  libmesh_assert(base_fe.get());
126  libmesh_assert(inf_elem);
127 
128  // checks for consistency of requested calculations,
129  // adds further quantities as needed.
130  this->determine_calculations();
131 
132  if (pts == nullptr)
133  {
134  libmesh_assert(base_fe->qrule);
135  libmesh_assert_equal_to (base_fe->qrule, base_qrule.get());
136  libmesh_assert(radial_qrule.get());
137 
138  bool init_shape_functions_required = false;
139 
140  // init the radial data fields only when the radial order changes
141  if (current_fe_type.radial_order != fe_type.radial_order)
142  {
143  current_fe_type.radial_order = fe_type.radial_order;
144 
145  // Watch out: this call to QBase->init() only works for
146  // current_fe_type = const! To allow variable Order,
147  // the init() of QBase has to be modified...
148  radial_qrule->init(EDGE2);
149 
150  // initialize the radial shape functions
151  this->init_radial_shape_functions(inf_elem);
152 
153  init_shape_functions_required=true;
154  }
155 
156 
157  bool update_base_elem_required=true;
158 
159  // update the type in accordance to the current cell
160  // and reinit if the cell type has changed or (as in
161  // the case of the hierarchics) the shape functions
162  // depend on the particular element and need a reinit
163  if ((Dim != 1) &&
164  ((this->get_type() != inf_elem->type()) ||
165  (base_fe->shapes_need_reinit())))
166  {
167  // store the new element type, update base_elem
168  // here. Through \p update_base_elem_required,
169  // remember whether it has to be updated (see below).
170  elem_type = inf_elem->type();
171  this->update_base_elem(inf_elem);
172  update_base_elem_required=false;
173 
174  // initialize the base quadrature rule for the new element
175  base_qrule->init(base_elem->type());
176  init_shape_functions_required=true;
177 
178  }
179 
180  // computing the reference-to-physical map and coordinates works
181  // only, if we have the current base_elem stored.
182  // This happens when fe_type is const,
183  // the inf_elem->type remains the same. Then we have to
184  // update the base elem _here_.
185  if (update_base_elem_required)
186  this->update_base_elem(inf_elem);
187 
188  if (calculate_phi_scaled || calculate_dphi_scaled || calculate_phi || calculate_dphi)
189  // initialize the shape functions in the base
190  base_fe->init_base_shape_functions(base_fe->qrule->get_points(),
191  base_elem.get());
192 
193  // compute the shape functions and map functions of base_fe
194  // before using them later in compute_shape_functions.
195  base_fe->_fe_map->compute_map (base_fe->dim, base_fe->qrule->get_weights(),
196  base_elem.get(), base_fe->calculate_d2phi);
197  base_fe->compute_shape_functions(base_elem.get(), base_fe->qrule->get_points());
198 
199  // when either the radial or base part change,
200  // we have to init the whole fields
201  if (init_shape_functions_required)
202  this->init_shape_functions (radial_qrule->get_points(),
203  base_fe->qrule->get_points(),
204  inf_elem);
205 
206  // Compute the shape functions and the derivatives
207  // at all quadrature points.
208  this->compute_shape_functions (inf_elem,
209  base_fe->qrule->get_points(),
210  radial_qrule->get_points()
211  /* weights are computed inside the function*/
212  );
213  }
214 
215  else // if pts != nullptr
216  {
217  // update the elem_type
218  elem_type = inf_elem->type();
219 
220  // We'll assume that pts is a tensor product mesh of points.
221  // pts[i] = pts[ angular_index + n_angular_pts * radial_index]
222  // That will handle the pts.size()==1 case that we care about
223  // right now, and it will generalize a bit, and it won't break
224  // the assumptions elsewhere in InfFE.
225  std::vector<Point> radial_pts;
226  if (pts->size() > 0)
227  {
228  Real radius = (*pts)[0](Dim-1);
229  radial_pts.push_back(radius);
230  unsigned int n_radial_pts=1;
231  unsigned int n_angular_pts=1;
232  for (auto p : IntRange<std::size_t>(1, pts->size()))
233  {
234  radius = (*pts)[p](Dim-1);
235  // check for changes of radius: The max. allowed distance is somewhat arbitrary
236  // but the given value should not produce false positives...
237  if (std::abs(radial_pts[n_radial_pts-1](0) - radius) > 1e-4)
238  {
239  // it may change only every n_angular_pts:
240  if (p == (n_radial_pts)*n_angular_pts)
241  {
242  radial_pts.push_back(radius);
243  ++n_radial_pts;
244  }
245  else
246  {
247  libmesh_error_msg("We assumed that the "<<pts->size()
248  <<" points are of tensor-product type with "
249  <<n_radial_pts<<" radial points and "
250  <<n_angular_pts<< " angular points."<<std::endl
251  <<"But apparently point "<<p+1
252  <<" does not fit that scheme: Its radius is "
253  <<radius <<"but should have "
254  <<radial_pts[n_radial_pts*n_angular_pts-p]<<".");
255  //<<radial_pts[p-n_radial_pts*n_angular_pts]<<".");
256  }
257  }
258  // if we are still at the first radial segment,
259  // we consider another angular point
260  else if (n_radial_pts == 1)
261  {
262  ++n_angular_pts;
263  }
264  // if there was repetition but this does not, the assumed
265  // format does not work:
266  }
267  }
268  else
269  {
270  // I don't see any reason to call this function with no points.
271  libmesh_error_msg("Calling reinit() with an empty point list is prohibited.\n");
272  }
273 
274  const std::size_t radial_pts_size = radial_pts.size();
275  const std::size_t base_pts_size = pts->size() / radial_pts_size;
276  // If we're a tensor product we should have no remainder
277  libmesh_assert_equal_to
278  (base_pts_size * radial_pts_size, pts->size());
279 
280 
281  std::vector<Point> base_pts;
282  base_pts.reserve(base_pts_size);
283  for (std::size_t p=0; p != base_pts_size; ++p)
284  {
285  Point pt = (*pts)[p];
286  pt(Dim-1) = 0;
287  base_pts.push_back(pt);
288  }
289 
290  // init radial shapes
291  this->init_radial_shape_functions(inf_elem, &radial_pts);
292 
293  // update the base
294  this->update_base_elem(inf_elem);
295 
296  // the finite element on the ifem base
297  base_fe = FEBase::build(Dim-1, this->fe_type);
298 
299  // having a new base_fe, we need to redetermine the tasks...
300  this->determine_calculations();
301 
302  base_fe->reinit( base_elem.get(), &base_pts);
303 
304  this->init_shape_functions (radial_pts, base_pts, inf_elem);
305 
306  // finally compute the ifem shapes
307  if (weights != nullptr)
308  {
309  this->compute_shape_functions (inf_elem,base_pts,radial_pts);
310  }
311  else
312  {
313  this->compute_shape_functions (inf_elem, base_pts, radial_pts);
314  }
315 
316  }
317 
318  if (this->calculate_dual)
319  libmesh_not_implemented_msg("Dual shape support for infinite elements is "
320  "not currently implemented");
321 }
322 
323 
324 
325 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
327 {
328  this->calculations_started = true;
329 
330  // If the user forgot to request anything, but we're enabling
331  // deprecated backwards compatibility, then we'll be safe and
332  // calculate everything. If we haven't enable deprecated backwards
333  // compatibility then we'll scream and die.
334 #ifdef LIBMESH_ENABLE_DEPRECATED
335  if (!this->calculate_nothing &&
336  !this->calculate_phi && !this->calculate_dphi &&
337  !this->calculate_dphiref &&
338  !this->calculate_phi_scaled && !this->calculate_dphi_scaled &&
339  !this->calculate_xyz && !this->calculate_jxw &&
340  !this->calculate_map_scaled && !this->calculate_map &&
341 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
342  !this->calculate_d2phi &&
343 #endif
344  !this->calculate_curl_phi && !this->calculate_div_phi)
345  {
346  libmesh_deprecated();
347  this->calculate_phi = this->calculate_dphi = this->calculate_jxw = true;
348  this->calculate_dphiref = true;
349 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
350  this->calculate_d2phi = true;
351 #endif
352  this->calculate_phi_scaled = this->calculate_dphi_scaled = this->calculate_xyz = true;
353  if (FEInterface::field_type(fe_type.family) == TYPE_VECTOR)
354  {
355  this->calculate_curl_phi = true;
356  this->calculate_div_phi = true;
357  }
358  }
359 #else //LIBMESH_ENABLE_DEPRECATED
360  // ANSI C does not allow to embed the preprocessor-statement into the assert, so we
361  // make two statements, just different by 'calculate_d2phi'.
362 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
363  libmesh_assert (this->calculate_nothing || this->calculate_d2phi ||
364  this->calculate_phi || this->calculate_dphi ||
365  this->calculate_dphiref ||
366  this->calculate_phi_scaled || this->calculate_dphi_scaled ||
367  this->calculate_xyz || this->calculate_jxw ||
368  this->calculate_map_scaled || this->calculate_map ||
369  this->calculate_curl_phi || this->calculate_div_phi);
370 #else
371  libmesh_assert (this->calculate_nothing ||
372  this->calculate_phi || this->calculate_dphi ||
373  this->calculate_dphiref ||
374  this->calculate_phi_scaled || this->calculate_dphi_scaled ||
375  this->calculate_xyz || this->calculate_jxw ||
376  this->calculate_map_scaled || this->calculate_map ||
377  this->calculate_curl_phi || this->calculate_div_phi);
378 #endif
379 #endif // LIBMESH_ENABLE_DEPRECATED
380 
381  // set further terms necessary to do the requested task
382  if (calculate_jxw)
383  this->calculate_map = true;
384  if (this->calculate_dphi)
385  this->calculate_map = true;
386  if (this->calculate_dphi_scaled)
387  this->calculate_map_scaled = true;
388  if (calculate_xyz && !calculate_map)
389  // if Cartesian positions were requested but the calculation of map
390  // was not triggered, we'll opt for the 'scaled' variant.
391  this->calculate_map_scaled = true;
392  base_fe->calculate_phi = this->calculate_phi || this->calculate_phi_scaled
393  || this->calculate_dphi || this->calculate_dphi_scaled;
394  base_fe->calculate_dphi = this->calculate_dphi || this->calculate_dphi_scaled;
395  if (this->calculate_map || this->calculate_map_scaled
396  || this->calculate_dphiref)
397  {
398  base_fe->calculate_dphiref = true;
399  base_fe->get_xyz(); // trigger base_fe->fe_map to 'calculate_xyz'
400  base_fe->get_JxW(); // trigger base_fe->fe_map to 'calculate_dxyz'
401  }
402  base_fe->determine_calculations();
403 
404 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
405  if (this->calculate_d2phi)
406  libmesh_not_implemented();
407 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES
408 }
409 
410 
411 
412 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
413 void
415 init_radial_shape_functions(const Elem * libmesh_dbg_var(inf_elem),
416  const std::vector<Point> * radial_pts)
417 {
418  libmesh_assert(radial_qrule.get() || radial_pts);
419  libmesh_assert(inf_elem);
420 
421  // Start logging the radial shape function initialization
422  LOG_SCOPE("init_radial_shape_functions()", "InfFE");
423 
424  // initialize most of the things related to physical approximation
425  const Order radial_approx_order = fe_type.radial_order;
426  const unsigned int n_radial_approx_shape_functions =
427  InfFERadial::n_dofs(radial_approx_order);
428 
429  const std::size_t n_radial_qp =
430  radial_pts ? radial_pts->size() : radial_qrule->n_points();
431  const std::vector<Point> & radial_qp =
432  radial_pts ? *radial_pts : radial_qrule->get_points();
433 
434  // the radial polynomials (eval)
435  if (calculate_phi || calculate_dphi || calculate_phi_scaled || calculate_dphi_scaled)
436  {
437  mode.resize (n_radial_approx_shape_functions);
438  for (unsigned int i=0; i<n_radial_approx_shape_functions; ++i)
439  mode[i].resize (n_radial_qp);
440 
441  // evaluate the mode shapes in radial direction at radial quadrature points
442  for (unsigned int i=0; i<n_radial_approx_shape_functions; ++i)
443  for (std::size_t p=0; p<n_radial_qp; ++p)
444  mode[i][p] = InfFE<Dim,T_radial,T_map>::eval (radial_qp[p](0), radial_approx_order, i);
445  }
446 
447  if (calculate_dphi || calculate_dphi_scaled)
448  {
449  dmodedv.resize (n_radial_approx_shape_functions);
450  for (unsigned int i=0; i<n_radial_approx_shape_functions; ++i)
451  dmodedv[i].resize (n_radial_qp);
452 
453  // evaluate the mode shapes in radial direction at radial quadrature points
454  for (unsigned int i=0; i<n_radial_approx_shape_functions; ++i)
455  for (std::size_t p=0; p<n_radial_qp; ++p)
456  dmodedv[i][p] = InfFE<Dim,T_radial,T_map>::eval_deriv (radial_qp[p](0), radial_approx_order, i);
457  }
458 
459  // the (1-v)/2 weight.
460  if (calculate_phi || calculate_phi_scaled || calculate_dphi || calculate_dphi_scaled)
461  {
462  som.resize (n_radial_qp);
463  // compute scalar values at radial quadrature points
464  for (std::size_t p=0; p<n_radial_qp; ++p)
465  som[p] = InfFERadial::decay (Dim, radial_qp[p](0));
466  }
467  if (calculate_dphi || calculate_dphi_scaled)
468  {
469  dsomdv.resize (n_radial_qp);
470  // compute scalar values at radial quadrature points
471  for (std::size_t p=0; p<n_radial_qp; ++p)
472  dsomdv[p] = InfFERadial::decay_deriv (Dim, radial_qp[p](0));
473  }
474 }
475 
476 
477 
478 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
479 void InfFE<Dim,T_radial,T_map>::init_shape_functions(const std::vector<Point> & radial_qp,
480  const std::vector<Point> & base_qp,
481  const Elem * inf_elem)
482 {
483  libmesh_assert(inf_elem);
484 
485  // Start logging the radial shape function initialization
486  LOG_SCOPE("init_shape_functions()", "InfFE");
487 
488  // fast access to some const ints for the radial data
489  const unsigned int n_radial_approx_sf = InfFERadial::n_dofs(fe_type.radial_order);
490  const std::size_t n_radial_qp = radial_qp.size();
491 #ifdef DEBUG
492  if (calculate_phi || calculate_dphi || calculate_phi_scaled || calculate_dphi_scaled)
493  libmesh_assert_equal_to(n_radial_approx_sf, mode.size());
494  if (calculate_phi || calculate_dphi || calculate_dphi_scaled)
495  libmesh_assert_equal_to(som.size(), n_radial_qp);
496 #endif
497 
498 
499  // initialize most of the quantities related to mapping
500 
501  // The element type and order to use in the base map
502  //const Order base_mapping_order = base_elem->default_order();
503 
504  // the number of base shape functions used to construct the map
505  // (Lagrange shape functions are used for mapping in the base)
506  //unsigned int n_base_mapping_shape_functions =
507  // InfFEBase::n_base_mapping_sf(*base_elem,
508  // base_mapping_order);
509 
510  // initialize most of the things related to physical approximation
511  unsigned int n_base_approx_shape_functions;
512  if (Dim > 1)
513  n_base_approx_shape_functions = base_fe->n_shape_functions();
514  else
515  n_base_approx_shape_functions = 1;
516 
517 
518  // update class member field
519  _n_total_approx_sf =
520  n_radial_approx_sf * n_base_approx_shape_functions;
521 
522 
523  // The number of the base quadrature points.
524  const unsigned int n_base_qp = cast_int<unsigned int>(base_qp.size());
525 
526  // The total number of quadrature points.
527  _n_total_qp = n_radial_qp * n_base_qp;
528 
529 
530  // initialize the node and shape numbering maps
531  {
532  // similar for the shapes: the i-th entry stores
533  // the associated base/radial shape number
534  _radial_shape_index.resize(_n_total_approx_sf);
535  _base_shape_index.resize(_n_total_approx_sf);
536 
537  // fill the shape index map
538  for (unsigned int n=0; n<_n_total_approx_sf; ++n)
539  {
540  compute_shape_indices (this->fe_type,
541  inf_elem,
542  n,
543  _base_shape_index[n],
544  _radial_shape_index[n]);
545  libmesh_assert_less (_base_shape_index[n], n_base_approx_shape_functions);
546  libmesh_assert_less (_radial_shape_index[n], n_radial_approx_sf);
547  }
548  }
549 
550  // resize the base data fields
551  //dist.resize(n_base_mapping_shape_functions);
552 
553  // resize the total data fields
554 
555  // the phase term varies with xi, eta and zeta(v): store it for _all_ qp
556  //
557  // when computing the phase, we need the base approximations
558  // therefore, initialize the phase here, but evaluate it
559  // in compute_shape_functions().
560  //
561  // the weight, though, is only needed at the radial quadrature points, n_radial_qp.
562  // but for a uniform interface to the protected data fields
563  // the weight data field (which are accessible from the outside) are expanded to _n_total_qp.
564  if (calculate_phi || calculate_dphi)
565  weight.resize (_n_total_qp);
566  if (calculate_phi_scaled || calculate_dphi_scaled)
567  weightxr_sq.resize (_n_total_qp);
568  if (calculate_dphi || calculate_dphi_scaled)
569  dweightdv.resize (n_radial_qp);
570  if (calculate_dphi)
571  dweight.resize (_n_total_qp);
572  if (calculate_dphi_scaled)
573  dweightxr_sq.resize(_n_total_qp);
574 
575  if (calculate_dphi || calculate_dphi_scaled)
576  dphase.resize (_n_total_qp);
577 
578  // this vector contains the integration weights for the combined quadrature rule
579  // if no quadrature rules are given, use only ones.
580  _total_qrule_weights.resize(_n_total_qp,1);
581 
582  // InfFE's data fields phi, dphi, dphidx, phi_map etc hold the _total_
583  // shape and mapping functions, respectively
584  {
585  if (calculate_map_scaled)
586  JxWxdecay.resize(_n_total_qp);
587  if (calculate_jxw)
588  JxW.resize(_n_total_qp);
589  if (calculate_map_scaled || calculate_map)
590  {
591  xyz.resize(_n_total_qp);
592  dxidx_map_scaled.resize(_n_total_qp);
593  dxidy_map_scaled.resize(_n_total_qp);
594  dxidz_map_scaled.resize(_n_total_qp);
595  detadx_map_scaled.resize(_n_total_qp);
596  detady_map_scaled.resize(_n_total_qp);
597  detadz_map_scaled.resize(_n_total_qp);
598  dzetadx_map_scaled.resize(_n_total_qp);
599  dzetady_map_scaled.resize(_n_total_qp);
600  dzetadz_map_scaled.resize(_n_total_qp);
601  }
602  if (calculate_map)
603  {
604  dxidx_map.resize(_n_total_qp);
605  dxidy_map.resize(_n_total_qp);
606  dxidz_map.resize(_n_total_qp);
607  detadx_map.resize(_n_total_qp);
608  detady_map.resize(_n_total_qp);
609  detadz_map.resize(_n_total_qp);
610  dzetadx_map.resize(_n_total_qp);
611  dzetady_map.resize(_n_total_qp);
612  dzetadz_map.resize(_n_total_qp);
613  }
614  if (calculate_phi)
615  phi.resize (_n_total_approx_sf);
616  if (calculate_phi_scaled)
617  phixr.resize (_n_total_approx_sf);
618  if (calculate_dphi)
619  {
620  dphi.resize (_n_total_approx_sf);
621  dphidx.resize (_n_total_approx_sf);
622  dphidy.resize (_n_total_approx_sf);
623  dphidz.resize (_n_total_approx_sf);
624  }
625 
626  if (calculate_dphi_scaled)
627  {
628  dphixr.resize (_n_total_approx_sf);
629  dphixr_sq.resize(_n_total_approx_sf);
630  }
631 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
632 
633  if (calculate_d2phi)
634  {
635  libmesh_not_implemented();
636  d2phi.resize (_n_total_approx_sf);
637  d2phidx2.resize (_n_total_approx_sf);
638  d2phidxdy.resize (_n_total_approx_sf);
639  d2phidxdz.resize (_n_total_approx_sf);
640  d2phidy2.resize (_n_total_approx_sf);
641  d2phidydz.resize (_n_total_approx_sf);
642  d2phidz2.resize (_n_total_approx_sf);
643  d2phidxi2.resize (_n_total_approx_sf);
644 
645  if (Dim > 1)
646  {
647  d2phidxideta.resize(_n_total_approx_sf);
648  d2phideta2.resize(_n_total_approx_sf);
649  }
650 
651  if (Dim > 2)
652  {
653  d2phidetadzeta.resize(_n_total_approx_sf);
654  d2phidxidzeta.resize(_n_total_approx_sf);
655  d2phidzeta2.resize(_n_total_approx_sf);
656  }
657  }
658 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
659 
660  if (calculate_dphi || calculate_dphi_scaled)
661  {
662  dphidxi.resize (_n_total_approx_sf);
663 
664  if (Dim > 1)
665  dphideta.resize(_n_total_approx_sf);
666 
667  if (Dim == 3)
668  dphidzeta.resize(_n_total_approx_sf);
669  }
670 
671  }
672 
673  // collect all the for loops, where inner vectors are
674  // resized to the appropriate number of quadrature points
675  {
676  if (calculate_phi)
677  for (unsigned int i=0; i<_n_total_approx_sf; ++i)
678  phi[i].resize (_n_total_qp);
679 
680  if (calculate_dphi)
681  for (unsigned int i=0; i<_n_total_approx_sf; ++i)
682  {
683  dphi[i].resize (_n_total_qp);
684  dphidx[i].resize (_n_total_qp);
685  dphidy[i].resize (_n_total_qp);
686  dphidz[i].resize (_n_total_qp);
687  }
688 
689  if (calculate_phi_scaled)
690  for (unsigned int i=0; i<_n_total_approx_sf; ++i)
691  {
692  phixr[i].resize (_n_total_qp);
693  }
694  if (calculate_dphi_scaled)
695  for (unsigned int i=0; i<_n_total_approx_sf; ++i)
696  {
697  dphixr[i].resize(_n_total_qp);
698  dphixr_sq[i].resize(_n_total_qp);
699  }
700 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
701  if (calculate_d2phi)
702  for (unsigned int i=0; i<_n_total_approx_sf; ++i)
703  {
704  d2phi[i].resize (_n_total_qp);
705  d2phidx2[i].resize (_n_total_qp);
706  d2phidxdy[i].resize (_n_total_qp);
707  d2phidxdz[i].resize (_n_total_qp);
708  d2phidy2[i].resize (_n_total_qp);
709  d2phidydz[i].resize (_n_total_qp);
710  d2phidy2[i].resize (_n_total_qp);
711  d2phidxi2[i].resize (_n_total_qp);
712 
713  if (Dim > 1)
714  {
715  d2phidxideta[i].resize (_n_total_qp);
716  d2phideta2[i].resize (_n_total_qp);
717  }
718  if (Dim > 2)
719  {
720  d2phidxidzeta[i].resize (_n_total_qp);
721  d2phidetadzeta[i].resize (_n_total_qp);
722  d2phidzeta2[i].resize (_n_total_qp);
723  }
724  }
725 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
726 
727  if (calculate_dphi || calculate_dphi_scaled)
728  for (unsigned int i=0; i<_n_total_approx_sf; ++i)
729  {
730  dphidxi[i].resize (_n_total_qp);
731 
732  if (Dim > 1)
733  dphideta[i].resize (_n_total_qp);
734 
735  if (Dim == 3)
736  dphidzeta[i].resize (_n_total_qp);
737 
738  }
739 
740  }
741  {
742  // (a) compute scalar values at _all_ quadrature points -- for uniform
743  // access from the outside to these fields
744  // (b) form a std::vector<Real> which contains the appropriate weights
745  // of the combined quadrature rule!
746  libmesh_assert_equal_to (radial_qp.size(), n_radial_qp);
747 
748  if (radial_qrule && base_qrule)
749  {
750  const std::vector<Real> & radial_qw = radial_qrule->get_weights();
751  const std::vector<Real> & base_qw = base_qrule->get_weights();
752  libmesh_assert_equal_to (radial_qw.size(), n_radial_qp);
753  libmesh_assert_equal_to (base_qw.size(), n_base_qp);
754 
755  for (unsigned int rp=0; rp<n_radial_qp; ++rp)
756  for (unsigned int bp=0; bp<n_base_qp; ++bp)
757  _total_qrule_weights[bp + rp*n_base_qp] = radial_qw[rp] * base_qw[bp];
758  }
759 
760 
761  for (unsigned int rp=0; rp<n_radial_qp; ++rp)
762  {
763  if (calculate_phi || calculate_dphi)
764  for (unsigned int bp=0; bp<n_base_qp; ++bp)
765  weight[bp + rp*n_base_qp] = InfFERadial::D(radial_qp[rp](0));
766 
767  if ( calculate_phi_scaled)
768  for (unsigned int bp=0; bp<n_base_qp; ++bp)
769  weightxr_sq[bp + rp*n_base_qp] = InfFERadial::Dxr_sq(radial_qp[rp](0));
770 
771  if ( calculate_dphi || calculate_dphi_scaled)
772  dweightdv[rp] = InfFERadial::D_deriv(radial_qp[rp](0));
773  }
774  }
775 }
776 
777 
778 
779 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
781  const std::vector<Point> & base_qp,
782  const std::vector<Point> & radial_qp
783  )
784 {
785  libmesh_assert(inf_elem);
786  // at least check whether the base element type is correct.
787  // otherwise this version of computing dist would give problems
788  libmesh_assert_equal_to (base_elem->type(),
789  InfFEBase::get_elem_type(inf_elem->type()));
790 
791  // Start logging the overall computation of shape functions
792  LOG_SCOPE("compute_shape_functions()", "InfFE");
793 
794  //const unsigned int n_radial_qp = cast_int<unsigned int>(som.size());
795  //const unsigned int n_base_qp = cast_int<unsigned int>(S_map[0].size());
796  const std::size_t n_radial_qp = radial_qp.size();
797  const unsigned int n_base_qp = base_qp.size();
798 
799  libmesh_assert_equal_to (_n_total_qp, n_radial_qp*n_base_qp);
800  libmesh_assert_equal_to (_n_total_qp, _total_qrule_weights.size());
801 #ifdef DEBUG
802  if (som.size() > 0)
803  libmesh_assert_equal_to(n_radial_qp, som.size());
804 
805  if (this->calculate_map || this->calculate_map_scaled)
806  {
807  // these vectors are needed later; initialize here already to have access to
808  // n_base_qp etc.
809  const std::vector<std::vector<Real>> & S_map = (base_fe->get_fe_map()).get_phi_map();
810  if (S_map[0].size() > 0)
811  libmesh_assert_equal_to(n_base_qp, S_map[0].size());
812  }
813  if (radial_qrule)
814  libmesh_assert_equal_to(n_radial_qp, radial_qrule->n_points());
815  if (base_qrule)
816  libmesh_assert_equal_to(n_base_qp, base_qrule->n_points());
817  libmesh_assert_equal_to(_n_total_qp % n_radial_qp, 0); // "Error in the structure of quadrature points!");
818 #endif
819 
820 
821  _n_total_approx_sf = InfFERadial::n_dofs(fe_type.radial_order) *
822  base_fe->n_shape_functions();
823 
824 
825 
826  const Point origin = inf_elem->origin();
827 
828  // Compute the shape function values (and derivatives)
829  // at the Quadrature points. Note that the actual values
830  // have already been computed via init_shape_functions
831 
832  unsigned int elem_dim = inf_elem->dim();
833  // Compute the value of the derivative shape function i at quadrature point p
834  switch (elem_dim)
835  {
836  case 1:
837  case 2:
838  {
839  libmesh_not_implemented();
840  break;
841  }
842  case 3:
843  {
844  std::vector<std::vector<Real>> S (0);
845  std::vector<std::vector<Real>> Ss(0);
846  std::vector<std::vector<Real>> St(0);
847 
848  std::vector<Real> base_dxidx (0);
849  std::vector<Real> base_dxidy (0);
850  std::vector<Real> base_dxidz (0);
851  std::vector<Real> base_detadx(0);
852  std::vector<Real> base_detady(0);
853  std::vector<Real> base_detadz(0);
854 
855  std::vector<Point> base_xyz (0);
856 
857  if (calculate_phi || calculate_phi_scaled ||
858  calculate_dphi || calculate_dphi_scaled)
859  S=base_fe->phi;
860 
861  // fast access to the approximation and mapping shapes of base_fe
862  if (calculate_map || calculate_map_scaled)
863  {
864  Ss = base_fe->dphidxi;
865  St = base_fe->dphideta;
866 
867  base_dxidx = base_fe->get_dxidx();
868  base_dxidy = base_fe->get_dxidy();
869  base_dxidz = base_fe->get_dxidz();
870  base_detadx = base_fe->get_detadx();
871  base_detady = base_fe->get_detady();
872  base_detadz = base_fe->get_detadz();
873 
874  base_xyz = base_fe->get_xyz();
875  }
876 
877  ElemType base_type= base_elem->type();
878 
879 #ifdef DEBUG
880  if (calculate_phi)
881  libmesh_assert_equal_to (phi.size(), _n_total_approx_sf);
882  if (calculate_dphi)
883  {
884  libmesh_assert_equal_to (dphidxi.size(), _n_total_approx_sf);
885  libmesh_assert_equal_to (dphideta.size(), _n_total_approx_sf);
886  libmesh_assert_equal_to (dphidzeta.size(), _n_total_approx_sf);
887  }
888 #endif
889 
890  unsigned int tp=0; // total qp
891  for (unsigned int rp=0; rp<n_radial_qp; ++rp) // over radial qps
892  for (unsigned int bp=0; bp<n_base_qp; ++bp) // over base qps
893 
894  { // First compute the map from base element quantities to physical space:
895 
896  // initialize them with invalid value to not use them
897  // without setting them to the correct value before.
898  Point unit_r(NAN);
899  RealGradient grad_a_scaled(NAN);
900  Real a(NAN);
901  Real r_norm(NAN);
902  if (calculate_map || calculate_map_scaled)
903  {
904  xyz[tp] = InfFEMap::map(elem_dim, inf_elem, Point(base_qp[bp](0),base_qp[bp](1),radial_qp[rp](0)));
905 
906  const Point r(xyz[tp]-origin);
907  a=(base_xyz[bp]-origin).norm();
908  r_norm = r.norm();
909 
910  // check that 'som' == a/r.
911 #ifndef NDEVEL
912  if (som.size())
913  libmesh_assert_less(std::abs(som[rp] -a/r_norm) , 1e-7);
914 #endif
915  unit_r=(r/r_norm);
916 
917  // They are used for computing the normal and do not correspond to the direction of eta and xi in this element:
918  // Due to the stretch of these axes in radial direction, they are deformed.
919  Point e_xi(base_dxidx[bp],
920  base_dxidy[bp],
921  base_dxidz[bp]);
922  Point e_eta(base_detadx[bp],
923  base_detady[bp],
924  base_detadz[bp]);
925 
926  const RealGradient normal=e_eta.cross(e_xi).unit();
927 
928  // grad a = a/r.norm() * grad_a_scaled
929  grad_a_scaled=unit_r - normal/(normal*unit_r);
930 
931  const Real dxi_er=base_dxidx[bp]* unit_r(0) + base_dxidy[bp] *unit_r(1) + base_dxidz[bp] *unit_r(2);
932  const Real deta_er=base_detadx[bp]*unit_r(0) + base_detady[bp]*unit_r(1) + base_detadz[bp]*unit_r(2);
933 
934  // in case of non-affine map, further terms need to be taken into account,
935  // involving \p e_eta and \p e_xi and thus recursive computation is needed
936  if (!base_elem->has_affine_map())
937  {
946  const unsigned int n_sf = base_elem->n_nodes();
947  RealGradient tmp(0.,0.,0.);
948  for (unsigned int i=0; i< n_sf; ++i)
949  {
950  RealGradient dL_da_i = (FE<2,LAGRANGE>::shape_deriv(base_type,
951  base_elem->default_order(),
952  i, 0, base_qp[bp]) * e_xi
953  +FE<2,LAGRANGE>::shape_deriv(base_type,
954  base_elem->default_order(),
955  i, 1, base_qp[bp]) * e_eta);
956 
957  tmp += (base_elem->node_ref(i) -origin).norm()* dL_da_i;
958 
959  }
960  libmesh_assert(tmp*unit_r < .95 ); // in a proper setup, tmp should have only a small radial component.
961  grad_a_scaled = ( tmp - (tmp*unit_r)*unit_r ) / ( 1. - tmp*unit_r);
962 
963  }
964 
965  // 'scale' = r/a
966  dxidx_map_scaled[tp] = (grad_a_scaled(0) - unit_r(0))*dxi_er +base_dxidx[bp];
967  dxidy_map_scaled[tp] = (grad_a_scaled(1) - unit_r(1))*dxi_er +base_dxidy[bp];
968  dxidz_map_scaled[tp] = (grad_a_scaled(2) - unit_r(2))*dxi_er +base_dxidz[bp];
969 
970  // 'scale' = r/a
971  detadx_map_scaled[tp] = (grad_a_scaled(0) - unit_r(0))*deta_er + base_detadx[bp];
972  detady_map_scaled[tp] = (grad_a_scaled(1) - unit_r(1))*deta_er + base_detady[bp];
973  detadz_map_scaled[tp] = (grad_a_scaled(2) - unit_r(2))*deta_er + base_detadz[bp];
974 
975  // 'scale' = (r/a)**2
976  dzetadx_map_scaled[tp] =-2./a*(grad_a_scaled(0) - unit_r(0));
977  dzetady_map_scaled[tp] =-2./a*(grad_a_scaled(1) - unit_r(1));
978  dzetadz_map_scaled[tp] =-2./a*(grad_a_scaled(2) - unit_r(2));
979 
980  }
981 
982  if (calculate_map)
983  {
984  dxidx_map[tp] = a/r_norm * dxidx_map_scaled[tp];
985  dxidy_map[tp] = a/r_norm * dxidy_map_scaled[tp];
986  dxidz_map[tp] = a/r_norm * dxidz_map_scaled[tp];
987 
988  detadx_map[tp] = a/r_norm * detadx_map_scaled[tp];
989  detady_map[tp] = a/r_norm * detady_map_scaled[tp];
990  detadz_map[tp] = a/r_norm * detadz_map_scaled[tp];
991 
992  // dzetadx = dzetadr*dr/dx - 2/r * grad_a
993  // = dzetadr*dr/dx - 2*a/r^2 * grad_a_scaled
994  dzetadx_map[tp] =-2.*a/(r_norm*r_norm)*(grad_a_scaled(0) - unit_r(0));
995  dzetady_map[tp] =-2.*a/(r_norm*r_norm)*(grad_a_scaled(1) - unit_r(1));
996  dzetadz_map[tp] =-2.*a/(r_norm*r_norm)*(grad_a_scaled(2) - unit_r(2));
997 
998  if (calculate_jxw)
999  {
1000  Real inv_jac = (dxidx_map[tp]*( detady_map[tp]*dzetadz_map[tp]- dzetady_map[tp]*detadz_map[tp]) +
1001  detadx_map[tp]*(dzetady_map[tp]* dxidz_map[tp]- dxidy_map[tp]*dzetadz_map[tp]) +
1002  dzetadx_map[tp]*( dxidy_map[tp]*detadz_map[tp]- detady_map[tp]* dxidz_map[tp]));
1003 
1004  if (inv_jac <= 1e-10)
1005  {
1006  libmesh_error_msg("ERROR: negative inverse Jacobian " \
1007  << inv_jac \
1008  << " at point " \
1009  << xyz[tp] \
1010  << " in element " \
1011  << inf_elem->id());
1012  }
1013 
1014 
1015  JxW[tp] = _total_qrule_weights[tp]/inv_jac;
1016  }
1017 
1018  }
1019  if (calculate_map_scaled)
1020  {
1021  Real inv_jacxR_pow4 = (dxidx_map_scaled[tp] *( detady_map_scaled[tp]*dzetadz_map_scaled[tp]
1022  - dzetady_map_scaled[tp]* detadz_map_scaled[tp]) +
1023  detadx_map_scaled[tp] *( dzetady_map_scaled[tp]* dxidz_map_scaled[tp]
1024  - dxidy_map_scaled[tp]*dzetadz_map_scaled[tp]) +
1025  dzetadx_map_scaled[tp]*( dxidy_map_scaled[tp]* detadz_map_scaled[tp]
1026  -detady_map_scaled[tp]* dxidz_map_scaled[tp]));
1027  if (inv_jacxR_pow4 <= 1e-7)
1028  {
1029  libmesh_error_msg("ERROR: negative weighted inverse Jacobian " \
1030  << inv_jacxR_pow4 \
1031  << " at point " \
1032  << xyz[tp] \
1033  << " in element " \
1034  << inf_elem->id());
1035  }
1036 
1037  JxWxdecay[tp] = _total_qrule_weights[tp]/inv_jacxR_pow4;
1038  }
1039 
1040  // phase term mu(r)=i*k*(r-a).
1041  // skip i*k: it is added separately during matrix assembly.
1042 
1043  if (calculate_dphi || calculate_dphi_scaled)
1044  dphase[tp] = unit_r - grad_a_scaled*a/r_norm;
1045 
1046  if (calculate_dphi)
1047  {
1048  dweight[tp](0) = dweightdv[rp] * dzetadx_map[tp];
1049  dweight[tp](1) = dweightdv[rp] * dzetady_map[tp];
1050  dweight[tp](2) = dweightdv[rp] * dzetadz_map[tp];
1051  }
1052  if (calculate_dphi_scaled)
1053  {
1054  dweightxr_sq[tp](0) = dweightdv[rp] * dzetadx_map_scaled[tp];
1055  dweightxr_sq[tp](1) = dweightdv[rp] * dzetady_map_scaled[tp];
1056  dweightxr_sq[tp](2) = dweightdv[rp] * dzetadz_map_scaled[tp];
1057  }
1058 
1059  if (calculate_phi || calculate_phi_scaled || calculate_dphi || calculate_dphi_scaled)
1060  // compute the shape-functions and derivative quantities:
1061  for (unsigned int i=0; i <_n_total_approx_sf ; ++i)
1062  {
1063  // let the index vectors take care of selecting the appropriate base/radial shape
1064  unsigned int bi = _base_shape_index [i];
1065  unsigned int ri = _radial_shape_index[i];
1066  if (calculate_phi)
1067  phi [i][tp] = S [bi][bp] * mode[ri][rp] * som[rp];
1068 
1069  if (calculate_phi_scaled)
1070  phixr [i][tp] = S [bi][bp] * mode[ri][rp];
1071 
1072  if (calculate_dphi || calculate_dphi_scaled)
1073  {
1074  dphidxi [i][tp] = Ss[bi][bp] * mode[ri][rp] * som[rp];
1075  dphideta [i][tp] = St[bi][bp] * mode[ri][rp] * som[rp];
1076  dphidzeta[i][tp] = S [bi][bp]
1077  * (dmodedv[ri][rp] * som[rp] + mode[ri][rp] * dsomdv[rp]);
1078  }
1079 
1080  if (calculate_dphi)
1081  {
1082 
1083  // dphi/dx = (dphi/dxi)*(dxi/dx) + (dphi/deta)*(deta/dx) + (dphi/dzeta)*(dzeta/dx);
1084  dphi[i][tp](0) =
1085  dphidx[i][tp] = (dphidxi[i][tp]*dxidx_map[tp] +
1086  dphideta[i][tp]*detadx_map[tp] +
1087  dphidzeta[i][tp]*dzetadx_map[tp]);
1088 
1089  // dphi/dy = (dphi/dxi)*(dxi/dy) + (dphi/deta)*(deta/dy) + (dphi/dzeta)*(dzeta/dy);
1090  dphi[i][tp](1) =
1091  dphidy[i][tp] = (dphidxi[i][tp]*dxidy_map[tp] +
1092  dphideta[i][tp]*detady_map[tp] +
1093  dphidzeta[i][tp]*dzetady_map[tp]);
1094 
1095  // dphi/dz = (dphi/dxi)*(dxi/dz) + (dphi/deta)*(deta/dz) + (dphi/dzeta)*(dzeta/dz);
1096  dphi[i][tp](2) =
1097  dphidz[i][tp] = (dphidxi[i][tp]*dxidz_map[tp] +
1098  dphideta[i][tp]*detadz_map[tp] +
1099  dphidzeta[i][tp]*dzetadz_map[tp]);
1100 
1101  }
1102  if (calculate_dphi_scaled)
1103  { // we don't distinguish between the different levels of scaling here...
1104 
1105  dphixr[i][tp](0)= (dphidxi[i][tp]*dxidx_map_scaled[tp] +
1106  dphideta[i][tp]*detadx_map_scaled[tp] +
1107  dphidzeta[i][tp]*dzetadx_map_scaled[tp]*som[rp]);
1108 
1109  dphixr[i][tp](1) = (dphidxi[i][tp]*dxidy_map_scaled[tp] +
1110  dphideta[i][tp]*detady_map_scaled[tp] +
1111  dphidzeta[i][tp]*dzetady_map_scaled[tp]*som[rp]);
1112 
1113  dphixr[i][tp](2) = (dphidxi[i][tp]*dxidz_map_scaled[tp] +
1114  dphideta[i][tp]*detadz_map_scaled[tp] +
1115  dphidzeta[i][tp]*dzetadz_map_scaled[tp]*som[rp]);
1116 
1117  const Real dphidxixr = Ss[bi][bp] * mode[ri][rp];
1118  const Real dphidetaxr= St[bi][bp] * mode[ri][rp];
1119 
1120  dphixr_sq[i][tp](0)= (dphidxixr*dxidx_map_scaled[tp] +
1121  dphidetaxr*detadx_map_scaled[tp] +
1122  dphidzeta[i][tp]*dzetadx_map_scaled[tp]);
1123 
1124  dphixr_sq[i][tp](1) = (dphidxixr*dxidy_map_scaled[tp] +
1125  dphidetaxr*detady_map_scaled[tp] +
1126  dphidzeta[i][tp]*dzetady_map_scaled[tp]);
1127 
1128  dphixr_sq[i][tp](2) = (dphidxixr*dxidz_map_scaled[tp] +
1129  dphidetaxr*detadz_map_scaled[tp] +
1130  dphidzeta[i][tp]*dzetadz_map_scaled[tp]);
1131  }
1132 
1133  }
1134  tp++;
1135  }
1136 
1137  break;
1138  }
1139 
1140  default:
1141  libmesh_error_msg("Unsupported dim = " << dim);
1142  }
1143 }
1144 
1145 
1146 
1147 template <unsigned int Dim, FEFamily T_radial, InfMapType T_map>
1149 {
1150  // We never call this.
1151  libmesh_not_implemented();
1152  return false;
1153 }
1154 
1155 } // namespace libMesh
1156 
1157 
1158 // Explicit instantiations
1159 #include "libmesh/inf_fe_instantiate_1D.h"
1160 #include "libmesh/inf_fe_instantiate_2D.h"
1161 #include "libmesh/inf_fe_instantiate_3D.h"
1162 
1163 
1164 
1165 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
Definition: fe_type.h:182
virtual void determine_calculations() override
Determine which values are to be calculated, for both the FE itself and for the FEMap.
Definition: inf_fe.C:326
ElemType
Defines an enum for geometric element types.
static ElemType get_elem_type(const ElemType type)
Order
defines an enum for polynomial orders.
Definition: enum_order.h:40
A specific instantiation of the FEBase class.
Definition: fe.h:42
static unsigned int n_dofs(const Order o_radial)
Definition: inf_fe.h:113
const Real radius
auto norm() const -> decltype(std::norm(T()))
Definition: type_vector.h:929
virtual Point origin() const
Definition: elem.h:1767
static Point map(const unsigned int dim, const Elem *inf_elem, const Point &reference_point)
Definition: inf_fe_map.C:40
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
static Real decay_deriv(const unsigned int dim, const Real)
Definition: inf_fe.h:1297
unsigned int dim
static Real Dxr_sq(const Real)
Definition: inf_fe.h:84
void init_shape_functions(const std::vector< Point > &radial_qp, const std::vector< Point > &base_qp, const Elem *inf_elem)
Initialize all the data fields like weight, mode, phi, dphidxi, dphideta, dphidzeta, etc.
Definition: inf_fe.C:479
This is the base class from which all geometric element types are derived.
Definition: elem.h:94
static FEFieldType field_type(const FEType &fe_type)
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
The libMesh namespace provides an interface to certain functionality in the library.
static std::unique_ptr< const Elem > build_elem(const Elem *inf_elem)
Build the base element of an infinite element.
virtual QuadratureType type() const =0
static Real D(const Real v)
Definition: inf_fe.h:82
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
Definition: type_vector.h:57
A specific instantiation of the FEBase class.
Definition: fe.h:127
friend class InfFE
Make all InfFE<Dim,T_radial,T_map> classes friends of each other, so that the protected eval() may be...
Definition: inf_fe.h:1258
dof_id_type id() const
Definition: dof_object.h:823
std::unique_ptr< FEBase > base_fe
Have a FE<Dim-1,T_base> handy for base approximation.
Definition: inf_fe.h:1211
dof_id_type weight(const MeshBase &mesh, const processor_id_type pid)
Definition: mesh_tools.C:436
static std::unique_ptr< FEGenericBase > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
libmesh_assert(ctx)
Order get_order() const
Definition: quadrature.h:218
static Real decay(const unsigned int dim, const Real v)
Definition: inf_fe.h:1271
unsigned int get_dim() const
Definition: quadrature.h:142
InfMapType inf_map
The coordinate mapping type of the infinite element.
Definition: fe_type.h:261
TypeVector< typename CompareTypes< T, T2 >::supertype > cross(const TypeVector< T2 > &v) const
Definition: type_vector.h:906
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
virtual bool shapes_need_reinit() const override
Definition: inf_fe.C:1148
void init_radial_shape_functions(const Elem *inf_elem, const std::vector< Point > *radial_pts=nullptr)
Some of the member data only depend on the radial part of the infinite element.
Definition: inf_fe.C:415
FEFamily radial_family
The type of approximation in radial direction.
Definition: fe_type.h:253
void compute_shape_functions(const Elem *inf_elem, const std::vector< Point > &base_qp, const std::vector< Point > &radial_qp)
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
Definition: inf_fe.C:780
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned short dim() const =0
static Real D_deriv(const Real v)
Definition: inf_fe.h:90
static std::unique_ptr< QBase > build(std::string_view name, const unsigned int dim, const Order order=INVALID_ORDER)
Builds a specific quadrature rule based on the name string.
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
Definition: inf_fe.C:121
virtual void attach_quadrature_rule(QBase *q) override
The use of quadrature rules with the InfFE class is somewhat different from the approach of the FE cl...
Definition: inf_fe.C:80
FEType fe_type
The finite element type for this object.
Definition: fe_abstract.h:709
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
The QBase class provides the basic functionality from which various quadrature rules can be derived...
Definition: quadrature.h:53
This class forms the foundation from which generic finite elements may be derived.
void update_base_elem(const Elem *inf_elem)
Updates the protected member base_elem to the appropriate base element for the given inf_elem...
Definition: inf_fe.C:110