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 //
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_COMMON_ARRAY_FACTORIES_IMPL_HPP
44 #define PANZER_COMMON_ARRAY_FACTORIES_IMPL_HPP
45 
46 #include "Teuchos_RCP.hpp"
47 #include "Phalanx_DataLayout_MDALayout.hpp"
48 #include "Phalanx_KokkosViewFactory.hpp"
49 
50 namespace panzer {
51 
52 // Implementation for intrepid container factory
53 template <typename Scalar,typename T0>
54 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
55 buildArray(const std::string & str,int d0) 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);
59 }
60 
61 template <typename Scalar,typename T0,typename T1>
62 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
63 buildArray(const std::string & str,int d0,int d1) const
64 {
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);
67 }
68 
69 template <typename Scalar,typename T0,typename T1,typename T2>
70 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
71 buildArray(const std::string & str,int d0,int d1,int d2) const
72 {
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);
75 }
76 
77 template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
78 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
79 buildArray(const std::string & str,int d0,int d1,int d2,int d3) const
80 {
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);
83 }
84 
85 template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
86 Kokkos::DynRankView<Scalar,PHX::Device> Intrepid2FieldContainerFactory::
87 buildArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
88 {
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);
91 }
92 
93 // Implementation for MDField array factory
94 template <typename Scalar,typename T0>
95 PHX::MDField<Scalar> MDFieldArrayFactory::
96 buildArray(const std::string & str,int d0) const
97 {
98  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
99 
100  PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0>(d0)));
101 
102  if(allocArray_)
103  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
104 
105  return field;
106 }
107 
108 template <typename Scalar,typename T0,typename T1>
109 PHX::MDField<Scalar> MDFieldArrayFactory::
110 buildArray(const std::string & str,int d0,int d1) const
111 {
112  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
113 
114  PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1 >(d0,d1)));
115 
116  if(allocArray_)
117  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
118 
119  return field;
120 }
121 
122 template <typename Scalar,typename T0,typename T1,typename T2>
123 PHX::MDField<Scalar> MDFieldArrayFactory::
124 buildArray(const std::string & str,int d0,int d1,int d2) const
125 {
126  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
127 
128  PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1,T2>(d0,d1,d2)));
129 
130  if(allocArray_)
131  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
132 
133  return field;
134 }
135 
136 template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
137 PHX::MDField<Scalar> MDFieldArrayFactory::
138 buildArray(const std::string & str,int d0,int d1,int d2,int d3) const
139 {
140  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
141 
142  PHX::MDField<Scalar> field = PHX::MDField<Scalar>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1,T2,T3>(d0,d1,d2,d3)));
143 
144  if(allocArray_)
145  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
146 
147  return field;
148 }
149 
150 template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
151 PHX::MDField<Scalar> MDFieldArrayFactory::
152 buildArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
153 {
154  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
155 
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)));
157 
158  if(allocArray_)
159  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
160 
161  return field;
162 }
163 
164 // Implementation for MDField array factory
165 template <typename Scalar,typename T0>
166 PHX::MDField<Scalar,T0> MDFieldArrayFactory::
167 buildStaticArray(const std::string & str,int d0) const
168 {
169  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
170 
171  PHX::MDField<Scalar,T0> field = PHX::MDField<Scalar,T0>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0>(d0)));
172 
173  if(allocArray_)
174  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
175 
176  return field;
177 }
178 
179 template <typename Scalar,typename T0,typename T1>
180 PHX::MDField<Scalar,T0,T1> MDFieldArrayFactory::
181 buildStaticArray(const std::string & str,int d0,int d1) const
182 {
183  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
184 
185  PHX::MDField<Scalar,T0,T1> field = PHX::MDField<Scalar,T0,T1>(prefix_+str,Teuchos::rcp(new PHX::MDALayout<T0,T1>(d0,d1)));
186 
187  if(allocArray_)
188  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
189 
190  return field;
191 }
192 
193 template <typename Scalar,typename T0,typename T1,typename T2>
194 PHX::MDField<Scalar,T0,T1,T2> MDFieldArrayFactory::
195 buildStaticArray(const std::string & str,int d0,int d1,int d2) const
196 {
197  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
198 
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)));
200 
201  if(allocArray_)
202  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
203 
204  return field;
205 }
206 
207 template <typename Scalar,typename T0,typename T1,typename T2,typename T3>
208 PHX::MDField<Scalar,T0,T1,T2,T3> MDFieldArrayFactory::
209 buildStaticArray(const std::string & str,int d0,int d1,int d2,int d3) const
210 {
211  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
212 
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)));
214 
215  if(allocArray_)
216  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
217 
218  return field;
219 }
220 
221 template <typename Scalar,typename T0,typename T1,typename T2,typename T3,typename T4>
222 PHX::MDField<Scalar,T0,T1,T2,T3,T4> MDFieldArrayFactory::
223 buildStaticArray(const std::string & str,int d0,int d1,int d2,int d3,int d4) const
224 {
225  typedef PHX::KokkosViewFactory<Scalar,PHX::Device> ViewFactory;
226 
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)));
228 
229  if(allocArray_)
230  field.setFieldData(ViewFactory::buildView(field.fieldTag(),ddims_));
231 
232  return field;
233 }
234 
235 } // end namespace panzer
236 
237 #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