Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_PointValues2_impl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #ifndef PANZER_POINT_VALUES2_IMPL_HPP
44 #define PANZER_POINT_VALUES2_IMPL_HPP
45 
46 #include "Intrepid2_CellTools.hpp"
47 
49 
50 // ***********************************************************
51 // * Evaluation and SetupArrays are NOT specialized
52 // ***********************************************************
53 
54 namespace panzer {
55 
56 
57 
58  template<typename Scalar,typename CoordinateArray>
59  struct CopyNodeCoords {
61  const CoordinateArray source_;
63  CopyNodeCoords(const CoordinateArray in_source,
65  : source_(in_source),target_(in_target) {}
66 
67  KOKKOS_INLINE_FUNCTION
68  void operator() (const size_type cell,const size_type node,const size_type dim) const {
69  target_(cell,node,dim) = source_(cell,node,dim);
70  }
71  };
72 
73  template<typename Scalar,typename CoordinateArray>
74  struct CopyPointCoords {
76  const CoordinateArray source_;
78  CopyPointCoords(const CoordinateArray in_source,
80  )
81  :source_(in_source),target_(in_target) {}
82 
83  KOKKOS_INLINE_FUNCTION
84  void operator() (const size_type pt,const size_type dim) const {
85  target_(pt,dim) = source_(pt,dim);
86  }
87  };
88 
89  template <typename Scalar>
92  {
93  MDFieldArrayFactory af(prefix_, ddims_, alloc_arrays_);
94 
95  point_rule = pr;
96 
97  int num_nodes = point_rule->topology->getNodeCount();
98  int num_cells = point_rule->workset_size;
99  int num_space_dim = point_rule->spatial_dimension;
100 
101  if (point_rule->isSide()) {
102  TEUCHOS_ASSERT(false); // not implemented!!!!
103  }
104 
105  int num_points = point_rule->num_points;
106 
107  coords_ref = af.template buildStaticArray<Scalar,IP,Dim>("coords_ref",num_points, num_space_dim);
108 
109  node_coordinates = af.template buildStaticArray<Scalar,Cell,NODE,Dim>("node_coordinates",num_cells, num_nodes, num_space_dim);
110 
111  jac = af.template buildStaticArray<Scalar,Cell,IP,Dim,Dim>("jac",num_cells, num_points, num_space_dim,num_space_dim);
112  jac_inv = af.template buildStaticArray<Scalar,Cell,IP,Dim,Dim>("jac_inv",num_cells, num_points, num_space_dim,num_space_dim);
113  jac_det = af.template buildStaticArray<Scalar,Cell,IP>("jac_det",num_cells, num_points);
114 
115  point_coords = af.template buildStaticArray<Scalar,Cell,IP,Dim>("point_coords",num_cells, num_points, num_space_dim);
116  }
117 
118  template <typename Scalar>
120  evaluateValues(const int in_num_cells)
121  {
122  if (point_rule->isSide()) {
123  TEUCHOS_ASSERT(false); // not implemented!!!!
124  }
125 
126  const int num_cells = in_num_cells < 0 ? (int) jac.extent(0) : in_num_cells;
127  const auto cell_range = std::pair<int,int>(0,num_cells);
128  auto s_jac = Kokkos::subview(jac.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
129  auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
130  auto s_jac_inv = Kokkos::subview(jac_inv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
131  auto s_node_coordinates = Kokkos::subview(node_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
132  auto s_point_coords = Kokkos::subview(point_coords.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
133  Intrepid2::CellTools<PHX::exec_space> cell_tools;
134 
135  cell_tools.setJacobian(s_jac, coords_ref.get_view(), s_node_coordinates, *(point_rule->topology));
136  cell_tools.setJacobianInv(s_jac_inv, s_jac);
137  cell_tools.setJacobianDet(s_jac_det, s_jac);
138 
139  // IP coordinates
140  cell_tools.mapToPhysicalFrame(s_point_coords, coords_ref.get_view(), s_node_coordinates, *(point_rule->topology));
141  }
142 
143  template <typename Scalar>
144  template <typename CoordinateArray>
146  copyNodeCoords(const CoordinateArray& in_node_coords)
147  {
148  // copy cell node coordinates
149  {
150  const size_type num_cells = in_node_coords.extent(0);
151  const size_type num_nodes = in_node_coords.extent(1);
152  const size_type num_dims = in_node_coords.extent(2);
153 
154  Kokkos::MDRangePolicy<PHX::Device::execution_space,Kokkos::Rank<3>> policy({0,0,0},{num_cells,num_nodes,num_dims});
155  Kokkos::parallel_for("PointValues2::copyNodeCoords",policy,panzer::CopyNodeCoords<Scalar,CoordinateArray>(in_node_coords,node_coordinates));
156  PHX::Device::execution_space().fence();
157  }
158  }
159 
160  template <typename Scalar>
161  template <typename CoordinateArray>
163  copyPointCoords(const CoordinateArray& in_point_coords)
164  {
165  // copy reference point values
166  {
167  const size_type num_points = in_point_coords.extent(0);
168  const size_type num_dims = in_point_coords.extent(1);
169 
170  Kokkos::MDRangePolicy<PHX::Device::execution_space,Kokkos::Rank<2>> policy({0,0},{num_points,num_dims});
171  Kokkos::parallel_for("PointValues2::copyPointCoords",policy,panzer::CopyPointCoords<Scalar,CoordinateArray>(in_point_coords,coords_ref));
172  PHX::Device::execution_space().fence();
173  }
174  }
175 
176 }
177 
178 #endif
typename panzer::PointValues2< Scalar >::size_type size_type
void evaluateValues(const CoordinateArray &node_coords, const PointArray &in_point_coords, const int in_num_cells=-1)
KOKKOS_INLINE_FUNCTION void operator()(const size_type cell, const size_type node, const size_type dim) const
PHX::MDField< Scalar, Cell, NODE, Dim > target_
PHX::MDField< Scalar, IP, Dim > target_
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
KOKKOS_INLINE_FUNCTION void operator()(const size_type pt, const size_type dim) const
typename panzer::PointValues2< 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)
CopyNodeCoords(const CoordinateArray in_source, PHX::MDField< Scalar, Cell, NODE, Dim > in_target)
CopyPointCoords(const CoordinateArray in_source, PHX::MDField< Scalar, IP, Dim > in_target)
#define TEUCHOS_ASSERT(assertion_test)