43 #ifndef PANZER_POINT_VALUES2_IMPL_HPP 
   44 #define PANZER_POINT_VALUES2_IMPL_HPP 
   46 #include "Intrepid2_CellTools.hpp" 
   56   template <
typename Scalar>
 
   64     int num_nodes = point_rule->
topology->getNodeCount();
 
   65     int num_cells = point_rule->workset_size;
 
   66     int num_space_dim = point_rule->spatial_dimension;
 
   68     if (point_rule->isSide()) {
 
   72     int num_points = point_rule->num_points;
 
   74     coords_ref = af.template buildStaticArray<Scalar,IP,Dim>(
"coords_ref",num_points, num_space_dim);
 
   76     node_coordinates = af.template buildStaticArray<Scalar,Cell,NODE,Dim>(
"node_coordinates",num_cells, num_nodes, num_space_dim);
 
   78     jac = af.template buildStaticArray<Scalar,Cell,IP,Dim,Dim>(
"jac",num_cells, num_points, num_space_dim,num_space_dim);
 
   79     jac_inv = af.template buildStaticArray<Scalar,Cell,IP,Dim,Dim>(
"jac_inv",num_cells, num_points, num_space_dim,num_space_dim);
 
   80     jac_det = af.template buildStaticArray<Scalar,Cell,IP>(
"jac_det",num_cells, num_points);
 
   82     point_coords = af.template buildStaticArray<Scalar,Cell,IP,Dim>(
"point_coords",num_cells, num_points, num_space_dim);
 
   85   template <
typename Scalar>
 
   89     if (point_rule->isSide()) {
 
   93     const int num_cells = in_num_cells < 0 ? (int) 
jac.extent(0) : in_num_cells;
 
   94     const auto cell_range = std::pair<int,int>(0,num_cells);
 
   95     auto s_jac = Kokkos::subview(
jac.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
 
   96     auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
 
   97     auto s_jac_inv = Kokkos::subview(jac_inv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
 
   98     auto s_node_coordinates = Kokkos::subview(node_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
 
   99     auto s_point_coords = Kokkos::subview(point_coords.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
 
  100     Intrepid2::CellTools<PHX::exec_space> cell_tools;
 
  102     cell_tools.setJacobian(s_jac, coords_ref.get_view(), s_node_coordinates, *(point_rule->topology));
 
  103     cell_tools.setJacobianInv(s_jac_inv, s_jac);
 
  104     cell_tools.setJacobianDet(s_jac_det, s_jac);
 
  107     cell_tools.mapToPhysicalFrame(s_point_coords, coords_ref.get_view(), s_node_coordinates, *(point_rule->topology));
 
  110   template <
typename Scalar>
 
  111   template <
typename CoordinateArray>
 
  117       size_type num_cells = in_node_coords.extent(0);
 
  118       size_type num_nodes = in_node_coords.extent(1);
 
  119       size_type num_dims = in_node_coords.extent(2);
 
  121       for (
size_type cell = 0; cell < num_cells;  ++cell)
 
  122   for (
size_type node = 0; node < num_nodes; ++node)
 
  123     for (
size_type dim = 0; dim < num_dims; ++dim)
 
  124       node_coordinates(cell,node,dim) = in_node_coords(cell,node,dim);
 
  128   template <
typename Scalar>
 
  129   template <
typename CoordinateArray>
 
  135       size_type num_points = in_point_coords.extent(0);
 
  136       size_type num_dims = in_point_coords.extent(1);
 
  138       for (
size_type point = 0; point < num_points; ++point)
 
  139         for (
size_type dim = 0; dim < num_dims; ++dim)
 
  140           coords_ref(point,dim) = in_point_coords(point,dim);
 
void evaluateValues(const CoordinateArray &node_coords, const PointArray &in_point_coords, const int in_num_cells=-1)
 
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays. 
 
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
 
void copyNodeCoords(const CoordinateArray &in_node_coords)
 
Teuchos::RCP< const shards::CellTopology > topology
 
void copyPointCoords(const CoordinateArray &in_point_coords)
 
#define TEUCHOS_ASSERT(assertion_test)