21 #include "libmesh/libmesh_config.h" 
   23 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 
   25 #include "libmesh/inf_fe.h" 
   26 #include "libmesh/quadrature_gauss.h" 
   27 #include "libmesh/elem.h" 
   28 #include "libmesh/libmesh_logging.h" 
   29 #include "libmesh/int_range.h" 
   30 #include "libmesh/auto_ptr.h" 
   38 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
   42   _n_total_approx_sf (0),
 
   52   current_fe_type (
FEType(fet.order,
 
   70 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
   77   const Order radial_int_order = static_cast<Order>(2 * (static_cast<unsigned int>(fe_type.radial_order.get_order()) + 1) +2);
 
   78   const unsigned int qrule_dim = q->
get_dim();
 
   84       base_fe->attach_quadrature_rule(base_qrule.get());
 
   88   radial_qrule = libmesh_make_unique<QGauss>(1, radial_int_order);
 
  100 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_base>
 
  111 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
  113                                        const std::vector<Point> * 
const pts,
 
  114                                        const std::vector<Real> * 
const weights)
 
  121   this->calculate_phi = this->calculate_dphi = this->calculate_dphiref = 
true;
 
  123   this->determine_calculations();
 
  124   base_fe->calculate_phi = base_fe->calculate_dphi = base_fe->calculate_dphiref = 
true;
 
  126   base_fe->determine_calculations();
 
  131       libmesh_assert_equal_to (base_fe->qrule, base_qrule.get());
 
  134       bool init_shape_functions_required = 
false;
 
  137       if (current_fe_type.radial_order != fe_type.radial_order)
 
  139           current_fe_type.radial_order = fe_type.radial_order;
 
  144           radial_qrule->init(
EDGE2);
 
  147           this->init_radial_shape_functions(inf_elem);
 
  149           init_shape_functions_required=
true;
 
  153       bool update_base_elem_required=
true;
 
  160           ((this->get_type() != inf_elem->
type())  ||
 
  161            (base_fe->shapes_need_reinit())))
 
  166           elem_type = inf_elem->
type();
 
  167           this->update_base_elem(inf_elem);
 
  168           update_base_elem_required=
false;
 
  171           base_qrule->init(base_elem->type());
 
  174           base_fe->init_base_shape_functions(base_fe->qrule->get_points(),
 
  179           base_fe->_fe_map->compute_map (base_fe->dim, base_fe->qrule->get_weights(),
 
  180                                          base_elem.get(), base_fe->calculate_d2phi);
 
  181           base_fe->compute_shape_functions(base_elem.get(), base_fe->qrule->get_points());
 
  183           init_shape_functions_required=
true;
 
  189       if (init_shape_functions_required)
 
  190         this->init_shape_functions (radial_qrule->get_points(),
 
  191                                     base_fe->qrule->get_points(),
 
  198       if (update_base_elem_required)
 
  199         this->update_base_elem(inf_elem);
 
  203       this->combine_base_radial (inf_elem);
 
  205       this->_fe_map->compute_map (this->
dim, _total_qrule_weights, inf_elem, this->calculate_d2phi);
 
  209       this->compute_shape_functions (inf_elem,base_fe->qrule->get_points());
 
  215       elem_type = inf_elem->
type();
 
  221       std::vector<Point> radial_pts;
 
  226           unsigned int n_radial_pts=1;
 
  227           unsigned int n_angular_pts=1;
 
  230               radius = (*pts)[p](Dim-1);
 
  233               if (
abs(radial_pts[p-n_radial_pts*n_angular_pts](0) - 
radius)< 1e-4)
 
  236                   if (p+1 == n_radial_pts*(n_angular_pts+1))
 
  240               else if (n_angular_pts == 1)
 
  246                   libmesh_error_msg(
"We assumed that the "<<pts->size()
 
  247                                   <<
" points are of tensor-product type with " 
  248                                   <<n_radial_pts<<
" radial points and " 
  249                                   <<n_angular_pts<< 
" angular points."<<std::endl
 
  250                                   <<
"But apparently point "<<p
 
  251                                   <<
" does not fit that scheme: Its radius is " 
  252                                   <<
radius <<
"but should have " 
  253                                   <<radial_pts[p-n_radial_pts*n_angular_pts]<<
".");
 
  257       const std::size_t radial_pts_size = radial_pts.size();
 
  258       const std::size_t base_pts_size = pts->size() / radial_pts_size;
 
  260       libmesh_assert_equal_to
 
  261         (base_pts_size * radial_pts_size, pts->size());
 
  264       std::vector<Point> base_pts;
 
  265       base_pts.reserve(base_pts_size);
 
  266       for (std::size_t p=0, ps=pts->size(); p != ps; p += radial_pts_size)
 
  268           Point pt = (*pts)[p];
 
  270           base_pts.push_back(pt);
 
  274       this->init_radial_shape_functions(inf_elem, &radial_pts);
 
  277       this->update_base_elem(inf_elem);
 
  282       base_fe->calculate_phi = base_fe->calculate_dphi = base_fe->calculate_dphiref = 
true;
 
  284       base_fe->determine_calculations();
 
  287       base_fe->init_base_shape_functions(base_pts,
 
  295           base_fe->_fe_map->compute_map (base_fe->dim, *weights,
 
  296                                          base_elem.get(), base_fe->calculate_d2phi);
 
  300           std::vector<Real> dummy_weights (pts->size(), 1.);
 
  301           base_fe->_fe_map->compute_map (base_fe->dim, dummy_weights,
 
  302                                          base_elem.get(), base_fe->calculate_d2phi);
 
  305       base_fe->compute_shape_functions(base_elem.get(), *pts);
 
  307       this->init_shape_functions (radial_pts, base_pts, inf_elem);
 
  310       this->combine_base_radial (inf_elem);
 
  313       if (weights != 
nullptr)
 
  315           this->_fe_map->compute_map (this->
dim, *weights, inf_elem, this->calculate_d2phi);
 
  319           std::vector<Real> dummy_weights (pts->size(), 1.);
 
  320           this->_fe_map->compute_map (this->
dim, dummy_weights, inf_elem, this->calculate_d2phi);
 
  324       this->compute_shape_functions (inf_elem,*pts);
 
  333 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
  337                             const std::vector<Point> * radial_pts)
 
  343   LOG_SCOPE(
"init_radial_shape_functions()", 
"InfFE");
 
  349   const unsigned int n_radial_mapping_shape_functions =
 
  353   const Order radial_approx_order = fe_type.radial_order;
 
  354   const unsigned int n_radial_approx_shape_functions =
 
  357   const std::size_t n_radial_qp =
 
  358     radial_pts ? radial_pts->size() : radial_qrule->n_points();
 
  359   const std::vector<Point> & radial_qp =
 
  360     radial_pts ? *radial_pts : radial_qrule->get_points();
 
  365   mode.resize      (n_radial_approx_shape_functions);
 
  366   dmodedv.resize   (n_radial_approx_shape_functions);
 
  369   som.resize       (n_radial_qp);
 
  370   dsomdv.resize    (n_radial_qp);
 
  373   radial_map.resize    (n_radial_mapping_shape_functions);
 
  374   dradialdv_map.resize (n_radial_mapping_shape_functions);
 
  377   for (
unsigned int i=0; i<n_radial_mapping_shape_functions; i++)
 
  379       radial_map[i].resize    (n_radial_qp);
 
  380       dradialdv_map[i].resize (n_radial_qp);
 
  384   for (
unsigned int i=0; i<n_radial_approx_shape_functions; i++)
 
  386       mode[i].resize    (n_radial_qp);
 
  387       dmodedv[i].resize (n_radial_qp);
 
  392   for (std::size_t p=0; p<n_radial_qp; p++)
 
  400   for (
unsigned int i=0; i<n_radial_approx_shape_functions; i++)
 
  401     for (std::size_t p=0; p<n_radial_qp; p++)
 
  409   for (
unsigned int i=0; i<n_radial_mapping_shape_functions; i++)
 
  410     for (std::size_t p=0; p<n_radial_qp; p++)
 
  412         radial_map[i][p] = 
InfFEMap::eval (radial_qp[p](0), radial_mapping_order, i);
 
  421 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
  423                                                      const std::vector<Point> & base_qp,
 
  424                                                      const Elem * inf_elem)
 
  429   LOG_SCOPE(
"init_shape_functions()", 
"InfFE");
 
  432   const unsigned int n_radial_mapping_sf = cast_int<unsigned int>(radial_map.size());
 
  433   const unsigned int n_radial_approx_sf  = cast_int<unsigned int>(mode.size());
 
  434   const unsigned int n_radial_qp         = cast_int<unsigned int>(som.size());
 
  440   const Order base_mapping_order = base_elem->default_order();
 
  444   unsigned int n_base_mapping_shape_functions =
 
  448   const unsigned int n_total_mapping_shape_functions =
 
  449     n_radial_mapping_sf * n_base_mapping_shape_functions;
 
  452   unsigned int n_base_approx_shape_functions;
 
  454     n_base_approx_shape_functions = base_fe->n_shape_functions();
 
  456     n_base_approx_shape_functions = 1;
 
  459   const unsigned int n_total_approx_shape_functions =
 
  460     n_radial_approx_sf * n_base_approx_shape_functions;
 
  463   _n_total_approx_sf = n_total_approx_shape_functions;
 
  467   const unsigned int n_base_qp = cast_int<unsigned int>(base_qp.size());
 
  470   const unsigned int n_total_qp = n_radial_qp * n_base_qp;
 
  474   _n_total_qp = n_total_qp;
 
  482     _radial_node_index.resize(n_total_mapping_shape_functions);
 
  483     _base_node_index.resize(n_total_mapping_shape_functions);
 
  487     _radial_shape_index.resize(n_total_approx_shape_functions);
 
  488     _base_shape_index.resize(n_total_approx_shape_functions);
 
  493     for (
unsigned int n=0; n<n_total_mapping_shape_functions; n++)
 
  495         compute_node_indices (inf_elem_type,
 
  498                               _radial_node_index[n]);
 
  499         libmesh_assert_less (_base_node_index[n], n_base_mapping_shape_functions);
 
  500         libmesh_assert_less (_radial_node_index[n], n_radial_mapping_sf);
 
  504     for (
unsigned int n=0; n<n_total_approx_shape_functions; n++)
 
  506         compute_shape_indices (this->fe_type,
 
  509                                _base_shape_index[n],
 
  510                                _radial_shape_index[n]);
 
  511         libmesh_assert_less (_base_shape_index[n], n_base_approx_shape_functions);
 
  512         libmesh_assert_less (_radial_shape_index[n], n_radial_approx_sf);
 
  517   dist.resize(n_base_mapping_shape_functions);
 
  530   weight.resize      (n_total_qp);
 
  531   dweightdv.resize   (n_total_qp);
 
  532   dweight.resize     (n_total_qp);
 
  534   dphase.resize      (n_total_qp);
 
  535   dphasedxi.resize   (n_total_qp);
 
  536   dphasedeta.resize  (n_total_qp);
 
  537   dphasedzeta.resize (n_total_qp);
 
  540   _total_qrule_weights.resize(n_total_qp);
 
  545     phi.resize     (n_total_approx_shape_functions);
 
  546     dphi.resize    (n_total_approx_shape_functions);
 
  547     dphidx.resize  (n_total_approx_shape_functions);
 
  548     dphidy.resize  (n_total_approx_shape_functions);
 
  549     dphidz.resize  (n_total_approx_shape_functions);
 
  550     dphidxi.resize (n_total_approx_shape_functions);
 
  551 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  552     libmesh_warning(
"Warning: Second derivatives for Infinite elements" 
  553                     << 
" are not yet implemented!" 
  556     d2phi.resize     (n_total_approx_shape_functions);
 
  557     d2phidx2.resize  (n_total_approx_shape_functions);
 
  558     d2phidxdy.resize (n_total_approx_shape_functions);
 
  559     d2phidxdz.resize (n_total_approx_shape_functions);
 
  560     d2phidy2.resize  (n_total_approx_shape_functions);
 
  561     d2phidydz.resize (n_total_approx_shape_functions);
 
  562     d2phidz2.resize  (n_total_approx_shape_functions);
 
  563     d2phidxi2.resize (n_total_approx_shape_functions);
 
  567         d2phidxideta.resize(n_total_approx_shape_functions);
 
  568         d2phideta2.resize(n_total_approx_shape_functions);
 
  573         d2phidetadzeta.resize(n_total_approx_shape_functions);
 
  574         d2phidxidzeta.resize(n_total_approx_shape_functions);
 
  575         d2phidzeta2.resize(n_total_approx_shape_functions);
 
  577 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  580       dphideta.resize(n_total_approx_shape_functions);
 
  583       dphidzeta.resize(n_total_approx_shape_functions);
 
  585     std::vector<std::vector<Real>> & phi_map = this->_fe_map->get_phi_map();
 
  586     std::vector<std::vector<Real>> & dphidxi_map = this->_fe_map->get_dphidxi_map();
 
  588     phi_map.resize(n_total_mapping_shape_functions);
 
  589     dphidxi_map.resize(n_total_mapping_shape_functions);
 
  590 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  591     std::vector<std::vector<Real>> & d2phidxi2_map = this->_fe_map->get_d2phidxi2_map();
 
  592     d2phidxi2_map.resize(n_total_mapping_shape_functions);
 
  596         std::vector<std::vector<Real>> & d2phidxideta_map = this->_fe_map->get_d2phidxideta_map();
 
  597         std::vector<std::vector<Real>> & d2phideta2_map = this->_fe_map->get_d2phideta2_map();
 
  598         d2phidxideta_map.resize(n_total_mapping_shape_functions);
 
  599         d2phideta2_map.resize(n_total_mapping_shape_functions);
 
  604         std::vector<std::vector<Real>> & d2phidxidzeta_map = this->_fe_map->get_d2phidxidzeta_map();
 
  605         std::vector<std::vector<Real>> & d2phidetadzeta_map = this->_fe_map->get_d2phidetadzeta_map();
 
  606         std::vector<std::vector<Real>> & d2phidzeta2_map = this->_fe_map->get_d2phidzeta2_map();
 
  607         d2phidxidzeta_map.resize(n_total_mapping_shape_functions);
 
  608         d2phidetadzeta_map.resize(n_total_mapping_shape_functions);
 
  609         d2phidzeta2_map.resize(n_total_mapping_shape_functions);
 
  611 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  615         std::vector<std::vector<Real>> & dphideta_map = this->_fe_map->get_dphideta_map();
 
  616         dphideta_map.resize(n_total_mapping_shape_functions);
 
  621         std::vector<std::vector<Real>> & dphidzeta_map = this->_fe_map->get_dphidzeta_map();
 
  622         dphidzeta_map.resize(n_total_mapping_shape_functions);
 
  629     for (
unsigned int i=0; i<n_total_approx_shape_functions; i++)
 
  631         phi[i].resize         (n_total_qp);
 
  632         dphi[i].resize        (n_total_qp);
 
  633         dphidx[i].resize      (n_total_qp);
 
  634         dphidy[i].resize      (n_total_qp);
 
  635         dphidz[i].resize      (n_total_qp);
 
  636         dphidxi[i].resize     (n_total_qp);
 
  637 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  638         d2phi[i].resize       (n_total_qp);
 
  639         d2phidx2[i].resize    (n_total_qp);
 
  640         d2phidxdy[i].resize   (n_total_qp);
 
  641         d2phidxdz[i].resize   (n_total_qp);
 
  642         d2phidy2[i].resize    (n_total_qp);
 
  643         d2phidydz[i].resize   (n_total_qp);
 
  644         d2phidy2[i].resize    (n_total_qp);
 
  645         d2phidxi2[i].resize   (n_total_qp);
 
  649             d2phidxideta[i].resize   (n_total_qp);
 
  650             d2phideta2[i].resize     (n_total_qp);
 
  654             d2phidxidzeta[i].resize  (n_total_qp);
 
  655             d2phidetadzeta[i].resize (n_total_qp);
 
  656             d2phidzeta2[i].resize    (n_total_qp);
 
  658 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  661           dphideta[i].resize  (n_total_qp);
 
  664           dphidzeta[i].resize (n_total_qp);
 
  668     for (
unsigned int i=0; i<n_total_mapping_shape_functions; i++)
 
  670         std::vector<std::vector<Real>> & phi_map = this->_fe_map->get_phi_map();
 
  671         std::vector<std::vector<Real>> & dphidxi_map = this->_fe_map->get_dphidxi_map();
 
  672         phi_map[i].resize         (n_total_qp);
 
  673         dphidxi_map[i].resize     (n_total_qp);
 
  674 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  675         std::vector<std::vector<Real>> & d2phidxi2_map = this->_fe_map->get_d2phidxi2_map();
 
  676         d2phidxi2_map[i].resize   (n_total_qp);
 
  679             std::vector<std::vector<Real>> & d2phidxideta_map = this->_fe_map->get_d2phidxideta_map();
 
  680             std::vector<std::vector<Real>> & d2phideta2_map = this->_fe_map->get_d2phideta2_map();
 
  681             d2phidxideta_map[i].resize   (n_total_qp);
 
  682             d2phideta2_map[i].resize     (n_total_qp);
 
  687             std::vector<std::vector<Real>> & d2phidxidzeta_map = this->_fe_map->get_d2phidxidzeta_map();
 
  688             std::vector<std::vector<Real>> & d2phidetadzeta_map = this->_fe_map->get_d2phidetadzeta_map();
 
  689             std::vector<std::vector<Real>> & d2phidzeta2_map = this->_fe_map->get_d2phidzeta2_map();
 
  690             d2phidxidzeta_map[i].resize  (n_total_qp);
 
  691             d2phidetadzeta_map[i].resize (n_total_qp);
 
  692             d2phidzeta2_map[i].resize    (n_total_qp);
 
  694 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 
  698             std::vector<std::vector<Real>> & dphideta_map = this->_fe_map->get_dphideta_map();
 
  699             dphideta_map[i].resize  (n_total_qp);
 
  704             std::vector<std::vector<Real>> & dphidzeta_map = this->_fe_map->get_dphidzeta_map();
 
  705             dphidzeta_map[i].resize (n_total_qp);
 
  717     libmesh_assert_equal_to (radial_qp.size(), n_radial_qp);
 
  719     if (radial_qrule && base_qrule)
 
  721         const std::vector<Real> & radial_qw = radial_qrule->get_weights();
 
  722         const std::vector<Real> & base_qw = base_qrule->get_weights();
 
  723         libmesh_assert_equal_to (radial_qw.size(), n_radial_qp);
 
  724         libmesh_assert_equal_to (base_qw.size(), n_base_qp);
 
  726         for (
unsigned int rp=0; rp<n_radial_qp; rp++)
 
  727           for (
unsigned int bp=0; bp<n_base_qp; bp++)
 
  731               _total_qrule_weights[bp + rp*n_base_qp] = radial_qw[rp] * base_qw[bp];
 
  736         for (
unsigned int rp=0; rp<n_radial_qp; rp++)
 
  737           for (
unsigned int bp=0; bp<n_base_qp; bp++)
 
  750 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
  756   libmesh_assert_equal_to (base_elem->type(),
 
  760   LOG_SCOPE(
"combine_base_radial()", 
"InfFE");
 
  763   std::fill (dphasedxi.begin(),   dphasedxi.end(),   0.);
 
  764   std::fill (dphasedeta.begin(),  dphasedeta.end(),  0.);
 
  765   std::fill (dphasedzeta.begin(), dphasedzeta.end(), 0.);
 
  768   const unsigned int n_base_mapping_sf = cast_int<unsigned int>(dist.size());
 
  772   for (
unsigned int n=0; n<n_base_mapping_sf; n++)
 
  773     dist[n] = 
Point(base_elem->point(n) - origin).
norm();
 
  781         libmesh_not_implemented();
 
  788         libmesh_not_implemented();
 
  796         const std::vector<std::vector<Real>> & S  = base_fe->phi;
 
  797         const std::vector<std::vector<Real>> & Ss = base_fe->dphidxi;
 
  798         const std::vector<std::vector<Real>> & St = base_fe->dphideta;
 
  799         const std::vector<std::vector<Real>> & S_map  = (base_fe->get_fe_map()).get_phi_map();
 
  800         const std::vector<std::vector<Real>> & Ss_map = (base_fe->get_fe_map()).get_dphidxi_map();
 
  801         const std::vector<std::vector<Real>> & St_map = (base_fe->get_fe_map()).get_dphideta_map();
 
  803         const unsigned int n_radial_qp = cast_int<unsigned int>(som.size());
 
  805           libmesh_assert_equal_to(n_radial_qp, radial_qrule->n_points());
 
  806         const unsigned int n_base_qp = cast_int<unsigned int>(S_map[0].size());
 
  808           libmesh_assert_equal_to(n_base_qp, base_qrule->n_points());
 
  810         const unsigned int n_total_mapping_sf =
 
  811           cast_int<unsigned int>(radial_map.size()) * n_base_mapping_sf;
 
  813         const unsigned int n_total_approx_sf =
 
  815           base_fe->n_shape_functions();
 
  821           for (
unsigned int rp=0; rp<n_radial_qp; rp++)  
 
  822             for (
unsigned int bp=0; bp<n_base_qp; bp++)  
 
  825                 for (
unsigned int i=0; i<n_base_mapping_sf; i++)
 
  827                     dphasedxi[tp]   += Ss_map[i][bp] * dist[i] * radial_map   [1][rp];
 
  828                     dphasedeta[tp]  += St_map[i][bp] * dist[i] * radial_map   [1][rp];
 
  829                     dphasedzeta[tp] += S_map [i][bp] * dist[i] * dradialdv_map[1][rp];
 
  837         libmesh_assert_equal_to (phi.size(), n_total_approx_sf);
 
  838         libmesh_assert_equal_to (dphidxi.size(), n_total_approx_sf);
 
  839         libmesh_assert_equal_to (dphideta.size(), n_total_approx_sf);
 
  840         libmesh_assert_equal_to (dphidzeta.size(), n_total_approx_sf);
 
  845         for (
unsigned int rp=0; rp<n_radial_qp; rp++)  
 
  846           for (
unsigned int bp=0; bp<n_base_qp; bp++)  
 
  847             for (
unsigned int ti=0; ti<n_total_approx_sf; ti++)  
 
  850                 const unsigned int bi = _base_shape_index  [ti];
 
  851                 const unsigned int ri = _radial_shape_index[ti];
 
  852                 phi      [ti][bp+rp*n_base_qp] = S [bi][bp] * mode[ri][rp] * som[rp];
 
  853                 dphidxi  [ti][bp+rp*n_base_qp] = Ss[bi][bp] * mode[ri][rp] * som[rp];
 
  854                 dphideta [ti][bp+rp*n_base_qp] = St[bi][bp] * mode[ri][rp] * som[rp];
 
  855                 dphidzeta[ti][bp+rp*n_base_qp] = S [bi][bp]
 
  856                   * (dmodedv[ri][rp] * som[rp] + mode[ri][rp] * dsomdv[rp]);
 
  859         std::vector<std::vector<Real>> & phi_map = this->_fe_map->get_phi_map();
 
  860         std::vector<std::vector<Real>> & dphidxi_map = this->_fe_map->get_dphidxi_map();
 
  861         std::vector<std::vector<Real>> & dphideta_map = this->_fe_map->get_dphideta_map();
 
  862         std::vector<std::vector<Real>> & dphidzeta_map = this->_fe_map->get_dphidzeta_map();
 
  864         libmesh_assert_equal_to (phi_map.size(), n_total_mapping_sf);
 
  865         libmesh_assert_equal_to (dphidxi_map.size(), n_total_mapping_sf);
 
  866         libmesh_assert_equal_to (dphideta_map.size(), n_total_mapping_sf);
 
  867         libmesh_assert_equal_to (dphidzeta_map.size(), n_total_mapping_sf);
 
  872         for (
unsigned int rp=0; rp<n_radial_qp; rp++) 
 
  873           for (
unsigned int bp=0; bp<n_base_qp; bp++) 
 
  874             for (
unsigned int ti=0; ti<n_total_mapping_sf; ti++) 
 
  877                 const unsigned int bi = _base_node_index  [ti];
 
  878                 const unsigned int ri = _radial_node_index[ti];
 
  879                 phi_map      [ti][bp+rp*n_base_qp] = S_map [bi][bp] * radial_map   [ri][rp];
 
  880                 dphidxi_map  [ti][bp+rp*n_base_qp] = Ss_map[bi][bp] * radial_map   [ri][rp];
 
  881                 dphideta_map [ti][bp+rp*n_base_qp] = St_map[bi][bp] * radial_map   [ri][rp];
 
  882                 dphidzeta_map[ti][bp+rp*n_base_qp] = S_map [bi][bp] * dradialdv_map[ri][rp];
 
  889       libmesh_error_msg(
"Unsupported Dim = " << Dim);
 
  898 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
  900                                                         const std::vector<Point> &)
 
  903   LOG_SCOPE(
"compute_shape_functions()", 
"InfFE");
 
  905   const unsigned int n_total_qp  = _n_total_qp;
 
  917         libmesh_not_implemented();
 
  923         libmesh_not_implemented();
 
  929         const std::vector<Real> & dxidx_map = this->_fe_map->get_dxidx();
 
  930         const std::vector<Real> & dxidy_map = this->_fe_map->get_dxidy();
 
  931         const std::vector<Real> & dxidz_map = this->_fe_map->get_dxidz();
 
  933         const std::vector<Real> & detadx_map = this->_fe_map->get_detadx();
 
  934         const std::vector<Real> & detady_map = this->_fe_map->get_detady();
 
  935         const std::vector<Real> & detadz_map = this->_fe_map->get_detadz();
 
  937         const std::vector<Real> & dzetadx_map = this->_fe_map->get_dzetadx();
 
  938         const std::vector<Real> & dzetady_map = this->_fe_map->get_dzetady();
 
  939         const std::vector<Real> & dzetadz_map = this->_fe_map->get_dzetadz();
 
  943           for (
unsigned int p=0; p<n_total_qp; p++)
 
  947                 dphidx[i][p] = (dphidxi[i][p]*dxidx_map[p] +
 
  948                                 dphideta[i][p]*detadx_map[p] +
 
  949                                 dphidzeta[i][p]*dzetadx_map[p]);
 
  953                 dphidy[i][p] = (dphidxi[i][p]*dxidy_map[p] +
 
  954                                 dphideta[i][p]*detady_map[p] +
 
  955                                 dphidzeta[i][p]*dzetady_map[p]);
 
  959                 dphidz[i][p] = (dphidxi[i][p]*dxidz_map[p] +
 
  960                                 dphideta[i][p]*detadz_map[p] +
 
  961                                 dphidzeta[i][p]*dzetadz_map[p]);
 
  966         for (
unsigned int p=0; p<n_total_qp; p++)
 
  969             dphase[p](0) = (dphasedxi[p]   * dxidx_map[p] +
 
  970                             dphasedeta[p]  * detadx_map[p] +
 
  971                             dphasedzeta[p] * dzetadx_map[p]);
 
  973             dphase[p](1) = (dphasedxi[p]   * dxidy_map[p] +
 
  974                             dphasedeta[p]  * detady_map[p] +
 
  975                             dphasedzeta[p] * dzetady_map[p]);
 
  977             dphase[p](2) = (dphasedxi[p]   * dxidz_map[p] +
 
  978                             dphasedeta[p]  * detadz_map[p] +
 
  979                             dphasedzeta[p] * dzetadz_map[p]);
 
  983             dweight[p](0) = dweightdv[p] * dzetadx_map[p];
 
  985             dweight[p](1) = dweightdv[p] * dzetady_map[p];
 
  987             dweight[p](2) = dweightdv[p] * dzetadz_map[p];
 
  994       libmesh_error_msg(
"Unsupported dim = " << 
dim);
 
 1000 template <
unsigned int Dim, FEFamily T_radial, InfMapType T_map>
 
 1010 #include "libmesh/inf_fe_instantiate_1D.h" 
 1011 #include "libmesh/inf_fe_instantiate_2D.h" 
 1012 #include "libmesh/inf_fe_instantiate_3D.h" 
 1016 #endif //ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS