43 #ifndef PANZER_COMMON_ARRAY_FACTORIES_IMPL_HPP
44 #define PANZER_COMMON_ARRAY_FACTORIES_IMPL_HPP
47 #include "Phalanx_DataLayout_MDALayout.hpp"
48 #include "Phalanx_KokkosViewFactory.hpp"
53 template <
typename Scalar,
typename T0>
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);
61 template <
typename Scalar,
typename T0,
typename T1>
65 static_assert(std::is_same<Scalar,double>::value,
"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
66 return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1);
69 template <
typename Scalar,
typename T0,
typename T1,
typename T2>
71 buildArray(
const std::string & str,
int d0,
int d1,
int d2)
const
73 static_assert(std::is_same<Scalar,double>::value,
"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
74 return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1,d2);
77 template <
typename Scalar,
typename T0,
typename T1,
typename T2,
typename T3>
79 buildArray(
const std::string & str,
int d0,
int d1,
int d2,
int d3)
const
81 static_assert(std::is_same<Scalar,double>::value,
"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
82 return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1,d2,d3);
85 template <
typename Scalar,
typename T0,
typename T1,
typename T2,
typename T3,
typename T4>
87 buildArray(
const std::string & str,
int d0,
int d1,
int d2,
int d3,
int d4)
const
89 static_assert(std::is_same<Scalar,double>::value,
"ERROR: CommonArryFactory for DynRankView only supports double scalar type!");
90 return Kokkos::DynRankView<Scalar,PHX::Device>(str,d0,d1,d2,d3,d4);
94 template <
typename Scalar,
typename T0>
98 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
103 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
108 template <
typename Scalar,
typename T0,
typename T1>
112 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
117 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
122 template <
typename Scalar,
typename T0,
typename T1,
typename T2>
124 buildArray(
const std::string & str,
int d0,
int d1,
int d2)
const
126 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
128 PHX::MDField<Scalar>
field = PHX::MDField<Scalar>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1,T2>(d0,d1,d2)));
131 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
136 template <
typename Scalar,
typename T0,
typename T1,
typename T2,
typename T3>
138 buildArray(
const std::string & str,
int d0,
int d1,
int d2,
int d3)
const
140 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
142 PHX::MDField<Scalar>
field = PHX::MDField<Scalar>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1,T2,T3>(d0,d1,d2,d3)));
145 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
150 template <
typename Scalar,
typename T0,
typename T1,
typename T2,
typename T3,
typename T4>
152 buildArray(
const std::string & str,
int d0,
int d1,
int d2,
int d3,
int d4)
const
154 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
156 PHX::MDField<Scalar>
field = PHX::MDField<Scalar>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1,T2,T3,T4>(d0,d1,d2,d3,d4)));
159 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
165 template <
typename Scalar,
typename T0>
169 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
174 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
179 template <
typename Scalar,
typename T0,
typename T1>
183 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
185 PHX::MDField<Scalar,T0,T1>
field = PHX::MDField<Scalar,T0,T1>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1>(d0,d1)));
188 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
193 template <
typename Scalar,
typename T0,
typename T1,
typename T2>
197 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
199 PHX::MDField<Scalar,T0,T1,T2>
field = PHX::MDField<Scalar,T0,T1,T2>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1,T2>(d0,d1,d2)));
202 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
207 template <
typename Scalar,
typename T0,
typename T1,
typename T2,
typename T3>
211 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
213 PHX::MDField<Scalar,T0,T1,T2,T3>
field = PHX::MDField<Scalar,T0,T1,T2,T3>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1,T2,T3>(d0,d1,d2,d3)));
216 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
221 template <
typename Scalar,
typename T0,
typename T1,
typename T2,
typename T3,
typename T4>
225 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
227 PHX::MDField<Scalar,T0,T1,T2,T3,T4>
field = PHX::MDField<Scalar,T0,T1,T2,T3,T4>(
prefix_+str,
Teuchos::rcp(
new PHX::MDALayout<T0,T1,T2,T3,T4>(d0,d1,d2,d3,d4)));
230 field.setFieldData(ViewFactory::buildView(field.fieldTag(),
ddims_));
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'll contribute, or in which we'll store, the result of computing this integral...
PHX::MDField< Scalar > buildArray(const std::string &str, int d0) const