Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_CommonArrayFactories_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_COMMON_ARRAY_FACTORIES_IMPL_HPP
12 #define PANZER_COMMON_ARRAY_FACTORIES_IMPL_HPP
13 
14 #include "Teuchos_RCP.hpp"
15 #include "Phalanx_DataLayout_MDALayout.hpp"
16 #include "Phalanx_KokkosViewFactory.hpp"
17 
18 namespace panzer {
19 
20 // Implementation for intrepid container factory
21 template <typename Scalar,typename T0>
22 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
23 buildArray(const std::string & str,int d0) const
24 {
25  static_assert(std::is_same<Scalar,double>::value,"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
26  return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0);
27 }
28 
29 template <typename Scalar,typename T0,typename T1>
30 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
31 buildArray(const std::string & str,int d0,int d1) const
32 {
33  static_assert(std::is_same<Scalar,double>::value,"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
34  return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1);
35 }
36 
37 template <typename Scalar,typename T0,typename T1,typename T2>
38 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
39 buildArray(const std::string & str,int d0,int d1,int d2) const
40 {
41  static_assert(std::is_same<Scalar,double>::value,"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
42  return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1,d2);
43 }
44 
45 template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
46 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
47 buildArray(const std::string & str,int d0,int d1,int d2,int d3) const
48 {
49  static_assert(std::is_same<Scalar,double>::value,"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
50  return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1,d2,d3);
51 }
52 
53 template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
54 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
55 buildArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
56 {
57  static_assert(std::is_same<Scalar,double>::value,"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
58  return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1,d2,d3,d4);
59 }
60 
61 // Implementation for MDField array factory
62 template <typename Scalar,typename T0>
64 buildArray(const std::string & str,int d0) const
65 {
67 
69 
70  if(allocArray_)
71  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
72 
73  return field;
74 }
75 
76 template <typename Scalar,typename T0,typename T1>
78 buildArray(const std::string & str,int d0,int d1) const
79 {
81 
83 
84  if(allocArray_)
85  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
86 
87  return field;
88 }
89 
90 template <typename Scalar,typename T0,typename T1,typename T2>
92 buildArray(const std::string & str,int d0,int d1,int d2) const
93 {
95 
97 
98  if(allocArray_)
99  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
100 
101  return field;
102 }
103 
104 template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
106 buildArray(const std::string & str,int d0,int d1,int d2,int d3) const
107 {
109 
111 
112  if(allocArray_)
113  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
114 
115  return field;
116 }
117 
118 template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
120 buildArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
121 {
123 
125 
126  if(allocArray_)
127  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
128 
129  return field;
130 }
131 
132 // Implementation for MDField array factory
133 template <typename Scalar,typename T0>
135 buildStaticArray(const std::string & str,int d0) const
136 {
138 
140 
141  if(allocArray_)
142  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
143 
144  return field;
145 }
146 
147 template <typename Scalar,typename T0,typename T1>
149 buildStaticArray(const std::string & str,int d0,int d1) const
150 {
152 
154 
155  if(allocArray_)
156  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
157 
158  return field;
159 }
160 
161 template <typename Scalar,typename T0,typename T1,typename T2>
163 buildStaticArray(const std::string & str,int d0,int d1,int d2) const
164 {
166 
168 
169  if(allocArray_)
170  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
171 
172  return field;
173 }
174 
175 template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
177 buildStaticArray(const std::string & str,int d0,int d1,int d2,int d3) const
178 {
180 
182 
183  if(allocArray_)
184  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
185 
186  return field;
187 }
188 
189 template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
191 buildStaticArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
192 {
194 
196 
197  if(allocArray_)
198  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
199 
200  return field;
201 }
202 
203 } // end namespace panzer
204 
205 #endif
std::vector< PHX::index_size_type > ddims_
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Kokkos::DynRankView< Scalar, PHX::Device > buildArray(const std::string &str, int d0) const
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
PHX::MDField< Scalar > buildArray(const std::string &str, int d0) const