Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_PointValues2_impl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_POINT_VALUES2_IMPL_HPP
12#define PANZER_POINT_VALUES2_IMPL_HPP
13
14#include "Intrepid2_CellTools.hpp"
15
17
18// ***********************************************************
19// * Evaluation and SetupArrays are NOT specialized
20// ***********************************************************
21
22namespace panzer {
23
24
25
26 template<typename Scalar,typename CoordinateArray>
29 const CoordinateArray source_;
30 PHX::MDField<Scalar,Cell,NODE,Dim> target_;
31 CopyNodeCoords(const CoordinateArray in_source,
32 PHX::MDField<Scalar,Cell,NODE,Dim> in_target)
33 : source_(in_source),target_(in_target) {}
34
35 KOKKOS_INLINE_FUNCTION
36 void operator() (const size_type cell,const size_type node,const size_type dim) const {
37 target_(cell,node,dim) = source_(cell,node,dim);
38 }
39 };
40
41 template<typename Scalar,typename CoordinateArray>
44 const CoordinateArray source_;
45 PHX::MDField<Scalar,IP,Dim> target_;
46 CopyPointCoords(const CoordinateArray in_source,
47 PHX::MDField<Scalar,IP,Dim> in_target
48 )
49 :source_(in_source),target_(in_target) {}
50
51 KOKKOS_INLINE_FUNCTION
52 void operator() (const size_type pt,const size_type dim) const {
53 target_(pt,dim) = source_(pt,dim);
54 }
55 };
56
57 template <typename Scalar>
59 setupArrays(const Teuchos::RCP<const PointRule> & pr)
60 {
61 MDFieldArrayFactory af(prefix_, ddims_, alloc_arrays_);
62
63 point_rule = pr;
64
65 int num_nodes = point_rule->topology->getNodeCount();
66 int num_cells = point_rule->workset_size;
67 int num_space_dim = point_rule->spatial_dimension;
68
69 if (point_rule->isSide()) {
70 TEUCHOS_ASSERT(false); // not implemented!!!!
71 }
72
73 int num_points = point_rule->num_points;
74
75 coords_ref = af.template buildStaticArray<Scalar,IP,Dim>("coords_ref",num_points, num_space_dim);
76
77 node_coordinates = af.template buildStaticArray<Scalar,Cell,NODE,Dim>("node_coordinates",num_cells, num_nodes, num_space_dim);
78
79 jac = af.template buildStaticArray<Scalar,Cell,IP,Dim,Dim>("jac",num_cells, num_points, num_space_dim,num_space_dim);
80 jac_inv = af.template buildStaticArray<Scalar,Cell,IP,Dim,Dim>("jac_inv",num_cells, num_points, num_space_dim,num_space_dim);
81 jac_det = af.template buildStaticArray<Scalar,Cell,IP>("jac_det",num_cells, num_points);
82
83 point_coords = af.template buildStaticArray<Scalar,Cell,IP,Dim>("point_coords",num_cells, num_points, num_space_dim);
84 }
85
86 template <typename Scalar>
88 evaluateValues(const int in_num_cells)
89 {
90 if (point_rule->isSide()) {
91 TEUCHOS_ASSERT(false); // not implemented!!!!
92 }
93
94 const int num_cells = in_num_cells < 0 ? (int) jac.extent(0) : in_num_cells;
95 const auto cell_range = std::pair<int,int>(0,num_cells);
96 auto s_jac = Kokkos::subview(jac.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
97 auto s_jac_det = Kokkos::subview(jac_det.get_view(),cell_range,Kokkos::ALL());
98 auto s_jac_inv = Kokkos::subview(jac_inv.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL());
99 auto s_node_coordinates = Kokkos::subview(node_coordinates.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
100 auto s_point_coords = Kokkos::subview(point_coords.get_view(),cell_range,Kokkos::ALL(),Kokkos::ALL());
101 Intrepid2::CellTools<PHX::exec_space> cell_tools;
102
103 cell_tools.setJacobian(s_jac, coords_ref.get_view(), s_node_coordinates, *(point_rule->topology));
104 cell_tools.setJacobianInv(s_jac_inv, s_jac);
105 cell_tools.setJacobianDet(s_jac_det, s_jac);
107 // IP coordinates
108 cell_tools.mapToPhysicalFrame(s_point_coords, coords_ref.get_view(), s_node_coordinates, *(point_rule->topology));
110
111 template <typename Scalar>
112 template <typename CoordinateArray>
114 copyNodeCoords(const CoordinateArray& in_node_coords)
115 {
116 // copy cell node coordinates
117 {
118 const size_type num_cells = in_node_coords.extent(0);
119 const size_type num_nodes = in_node_coords.extent(1);
120 const size_type num_dims = in_node_coords.extent(2);
121
122 Kokkos::MDRangePolicy<PHX::Device::execution_space,Kokkos::Rank<3>> policy({0,0,0},{num_cells,num_nodes,num_dims});
123 Kokkos::parallel_for("PointValues2::copyNodeCoords",policy,panzer::CopyNodeCoords<Scalar,CoordinateArray>(in_node_coords,node_coordinates));
124 PHX::Device::execution_space().fence();
125 }
126 }
127
128 template <typename Scalar>
129 template <typename CoordinateArray>
131 copyPointCoords(const CoordinateArray& in_point_coords)
132 {
133 // copy reference point values
134 {
135 const size_type num_points = in_point_coords.extent(0);
136 const size_type num_dims = in_point_coords.extent(1);
137
138 Kokkos::MDRangePolicy<PHX::Device::execution_space,Kokkos::Rank<2>> policy({0,0},{num_points,num_dims});
139 Kokkos::parallel_for("PointValues2::copyPointCoords",policy,panzer::CopyPointCoords<Scalar,CoordinateArray>(in_point_coords,coords_ref));
140 PHX::Device::execution_space().fence();
141 }
142 }
143
144}
145
146#endif
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
void copyPointCoords(const CoordinateArray &in_point_coords)
void evaluateValues(const CoordinateArray &node_coords, const PointArray &in_point_coords, const int in_num_cells=-1)
void copyNodeCoords(const CoordinateArray &in_node_coords)
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
KOKKOS_INLINE_FUNCTION void operator()(const size_type cell, const size_type node, const size_type dim) const
CopyNodeCoords(const CoordinateArray in_source, PHX::MDField< Scalar, Cell, NODE, Dim > in_target)
PHX::MDField< Scalar, Cell, NODE, Dim > target_
typename panzer::PointValues2< Scalar >::size_type size_type
CopyPointCoords(const CoordinateArray in_source, PHX::MDField< Scalar, IP, Dim > in_target)
typename panzer::PointValues2< Scalar >::size_type size_type
PHX::MDField< Scalar, IP, Dim > target_
KOKKOS_INLINE_FUNCTION void operator()(const size_type pt, const size_type dim) const