Panzer Version of the Day
Loading...
Searching...
No Matches
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
18namespace panzer {
19
20// Implementation for intrepid container factory
21template <typename Scalar,typename T0>
22Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
23buildArray(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
29template <typename Scalar,typename T0,typename T1>
30Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
31buildArray(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
37template <typename Scalar,typename T0,typename T1,typename T2>
38Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
39buildArray(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
45template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
46Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
47buildArray(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
53template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
54Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
55buildArray(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
62template <typename Scalar,typename T0>
63PHX::MDField<Scalar> MDFieldArrayFactory::
64buildArray(const std::string & str,int d0) const
65{
66 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
67
68 PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0>(d0)));
69
70 if(allocArray_)
71 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
72
73 return field;
74}
75
76template <typename Scalar,typename T0,typename T1>
77PHX::MDField<Scalar> MDFieldArrayFactory::
78buildArray(const std::string & str,int d0,int d1) const
79{
80 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
81
82 PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1 >(d0,d1)));
83
84 if(allocArray_)
85 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
86
87 return field;
88}
89
90template <typename Scalar,typename T0,typename T1,typename T2>
91PHX::MDField<Scalar> MDFieldArrayFactory::
92buildArray(const std::string & str,int d0,int d1,int d2) const
93{
94 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
95
96 PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1,T2>(d0,d1,d2)));
97
98 if(allocArray_)
99 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
100
101 return field;
102}
103
104template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
105PHX::MDField<Scalar> MDFieldArrayFactory::
106buildArray(const std::string & str,int d0,int d1,int d2,int d3) const
107{
108 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
109
110 PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1,T2,T3>(d0,d1,d2,d3)));
111
112 if(allocArray_)
113 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
114
115 return field;
116}
117
118template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
119PHX::MDField<Scalar> MDFieldArrayFactory::
120buildArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
121{
122 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
123
124 PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1,T2,T3,T4>(d0,d1,d2,d3,d4)));
125
126 if(allocArray_)
127 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
128
129 return field;
130}
131
132// Implementation for MDField array factory
133template <typename Scalar,typename T0>
134PHX::MDField<Scalar,T0> MDFieldArrayFactory::
135buildStaticArray(const std::string & str,int d0) const
136{
137 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
138
139 PHX::MDField<Scalar,T0> field = PHX::MDField<Scalar,T0>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0>(d0)));
140
141 if(allocArray_)
142 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
143
144 return field;
145}
146
147template <typename Scalar,typename T0,typename T1>
148PHX::MDField<Scalar,T0,T1> MDFieldArrayFactory::
149buildStaticArray(const std::string & str,int d0,int d1) const
150{
151 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
152
153 PHX::MDField<Scalar,T0,T1> field = PHX::MDField<Scalar,T0,T1>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1>(d0,d1)));
154
155 if(allocArray_)
156 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
157
158 return field;
159}
160
161template <typename Scalar,typename T0,typename T1,typename T2>
162PHX::MDField<Scalar,T0,T1,T2> MDFieldArrayFactory::
163buildStaticArray(const std::string & str,int d0,int d1,int d2) const
164{
165 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
166
167 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)));
168
169 if(allocArray_)
170 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
171
172 return field;
173}
174
175template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
176PHX::MDField<Scalar,T0,T1,T2,T3> MDFieldArrayFactory::
177buildStaticArray(const std::string & str,int d0,int d1,int d2,int d3) const
178{
179 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
180
181 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)));
182
183 if(allocArray_)
184 field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
185
186 return field;
187}
188
189template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
190PHX::MDField<Scalar,T0,T1,T2,T3,T4> MDFieldArrayFactory::
191buildStaticArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
192{
193 typedef PHX::KokkosViewFactory<Scalar,typename PHX::DevLayout<Scalar>::type,PHX::Device> ViewFactory;
194
195 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)));
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
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.
Kokkos::DynRankView< Scalar, PHX::Device > buildArray(const std::string &str, int d0) const
std::vector< PHX::index_size_type > ddims_
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
PHX::MDField< Scalar > buildArray(const std::string &str, int d0) const