49 #ifndef __INTREPID2_HGRAD_HEX_CN_FEMDEF_HPP__
50 #define __INTREPID2_HGRAD_HEX_CN_FEMDEF_HPP__
57 template<EOperator opType>
58 template<
typename OutputViewType,
59 typename inputViewType,
60 typename workViewType,
61 typename vinvViewType>
62 KOKKOS_INLINE_FUNCTION
64 Basis_HGRAD_HEX_Cn_FEM::Serial<opType>::
65 getValues( OutputViewType output,
66 const inputViewType input,
68 const vinvViewType vinv,
69 const ordinal_type operatorDn ) {
70 ordinal_type opDn = operatorDn;
72 const ordinal_type cardLine = vinv.extent(0);
73 const ordinal_type npts = input.extent(0);
75 typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
76 const auto input_x = Kokkos::subview(input, Kokkos::ALL(), range_type(0,1));
77 const auto input_y = Kokkos::subview(input, Kokkos::ALL(), range_type(1,2));
78 const auto input_z = Kokkos::subview(input, Kokkos::ALL(), range_type(2,3));
80 const ordinal_type dim_s = get_dimension_scalar(work);
81 auto ptr0 = work.data();
82 auto ptr1 = work.data()+cardLine*npts*dim_s;
83 auto ptr2 = work.data()+2*cardLine*npts*dim_s;
84 auto ptr3 = work.data()+3*cardLine*npts*dim_s;
86 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
87 auto vcprop = Kokkos::common_view_alloc_prop(work);
90 case OPERATOR_VALUE: {
91 viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
92 viewType output_x(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
93 viewType output_y(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
94 viewType output_z(Kokkos::view_wrap(ptr3, vcprop), cardLine, npts);
96 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
97 getValues(output_x, input_x, work_line, vinv);
99 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
100 getValues(output_y, input_y, work_line, vinv);
102 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
103 getValues(output_z, input_z, work_line, vinv);
106 ordinal_type idx = 0;
107 for (ordinal_type k=0;k<cardLine;++k)
108 for (ordinal_type j=0;j<cardLine;++j)
109 for (ordinal_type i=0;i<cardLine;++i,++idx)
110 for (ordinal_type l=0;l<npts;++l)
111 output.access(idx,l) = output_x.access(i,l)*output_y.access(j,l)*output_z.access(k,l);
125 opDn = getOperatorOrder(opType);
127 const ordinal_type dkcard = opDn + 1;
130 for (ordinal_type l1=0;l1<dkcard;++l1)
131 for (ordinal_type l0=0;l0<(l1+1);++l0) {
132 const ordinal_type mult_x = (opDn - l1);
133 const ordinal_type mult_y = l1 - l0;
134 const ordinal_type mult_z = l0;
142 viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
143 decltype(work_line) output_x, output_y, output_z;
146 output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts, 1);
147 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
148 getValues(output_x, input_x, work_line, vinv, mult_x);
150 output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
151 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
152 getValues(output_x, input_x, work_line, vinv);
156 output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts, 1);
157 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
158 getValues(output_y, input_y, work_line, vinv, mult_y);
160 output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
161 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
162 getValues(output_y, input_y, work_line, vinv);
166 output_z = viewType(Kokkos::view_wrap(ptr3, vcprop), cardLine, npts, 1);
167 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
168 getValues(output_z, input_z, work_line, vinv, mult_z);
170 output_z = viewType(Kokkos::view_wrap(ptr3, vcprop), cardLine, npts);
171 Impl::Basis_HGRAD_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
172 getValues(output_z, input_z, work_line, vinv);
176 ordinal_type idx = 0;
177 for (ordinal_type k=0;k<cardLine;++k)
178 for (ordinal_type j=0;j<cardLine;++j)
179 for (ordinal_type i=0;i<cardLine;++i,++idx)
180 for (ordinal_type l=0;l<npts;++l)
181 output.access(idx,l,d) = output_x.access(i,l,0)*output_y.access(j,l,0)*output_z.access(k,l,0);
188 INTREPID2_TEST_FOR_ABORT(
true ,
189 ">>> ERROR (Basis_HGRAD_HEX_Cn_FEM): Operator type not implemented");
195 template<
typename SpT, ordinal_type numPtsPerEval,
196 typename outputValueValueType,
class ...outputValueProperties,
197 typename inputPointValueType,
class ...inputPointProperties,
198 typename vinvValueType,
class ...vinvProperties>
200 Basis_HGRAD_HEX_Cn_FEM::
201 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
202 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
203 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
204 const EOperator operatorType ) {
205 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
206 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
207 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
208 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
211 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
212 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
213 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
214 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
216 typedef typename inputPointViewType::value_type inputPointType;
218 const ordinal_type cardinality = outputValues.extent(0);
219 const ordinal_type cardLine = std::cbrt(cardinality);
220 const ordinal_type workSize = 4*cardLine;
222 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
223 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
224 workViewType work(Kokkos::view_alloc(
"Basis_HGRAD_HEX_Cn_FEM::getValues::work", vcprop), workSize, inputPoints.extent(0));
226 switch (operatorType) {
227 case OPERATOR_VALUE: {
228 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
229 OPERATOR_VALUE,numPtsPerEval> FunctorType;
230 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
233 case OPERATOR_CURL: {
234 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
235 OPERATOR_CURL,numPtsPerEval> FunctorType;
236 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
250 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
251 OPERATOR_Dn,numPtsPerEval> FunctorType;
252 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
253 getOperatorOrder(operatorType)) );
257 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
258 ">>> ERROR (Basis_HGRAD_HEX_Cn_FEM): Operator type not implemented" );
266 template<
typename SpT,
typename OT,
typename PT>
269 const EPointType pointType ) {
275 this->vinv_ = Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>(
"Hgrad::HEX::Cn::vinv", cardLine, cardLine);
276 lineBasis.getVandermondeInverse(this->vinv_);
278 this->basisCardinality_ = cardLine*cardLine*cardLine;
279 this->basisDegree_ = order;
280 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<8> >() );
281 this->basisType_ = BASIS_FEM_FIAT;
282 this->basisCoordinates_ = COORDINATES_CARTESIAN;
283 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
288 const ordinal_type tagSize = 4;
289 const ordinal_type posScDim = 0;
290 const ordinal_type posScOrd = 1;
291 const ordinal_type posDfOrd = 2;
295 ordinal_type tags[maxCardLine*maxCardLine*maxCardLine][4];
297 const ordinal_type vert[2][2][2] = { { {0,1}, {3,2} },
300 const ordinal_type edge_x[2][2] = { {0, 4}, {2, 6} };
301 const ordinal_type edge_y[2][2] = { {3, 7}, {1, 5} };
302 const ordinal_type edge_z[2][2] = { {8,11}, {9,10} };
304 const ordinal_type face_yz[2] = {3, 1};
305 const ordinal_type face_xz[2] = {0, 2};
306 const ordinal_type face_xy[2] = {4, 5};
309 ordinal_type idx = 0;
310 for (
auto k=0;k<cardLine;++k) {
311 const auto tag_z = lineBasis.
getDofTag(k);
312 for (ordinal_type j=0;j<cardLine;++j) {
313 const auto tag_y = lineBasis.
getDofTag(j);
314 for (ordinal_type i=0;i<cardLine;++i,++idx) {
315 const auto tag_x = lineBasis.
getDofTag(i);
317 if (tag_x(0) == 0 && tag_y(0) == 0 && tag_z(0) == 0) {
320 tags[idx][1] = vert[tag_z(1)][tag_y(1)][tag_x(1)];
323 }
else if (tag_x(0) == 1 && tag_y(0) == 0 && tag_z(0) == 0) {
326 tags[idx][1] = edge_x[tag_y(1)][tag_z(1)];
327 tags[idx][2] = tag_x(2);
328 tags[idx][3] = tag_x(3);
329 }
else if (tag_x(0) == 0 && tag_y(0) == 1 && tag_z(0) == 0) {
332 tags[idx][1] = edge_y[tag_x(1)][tag_z(1)];
333 tags[idx][2] = tag_y(2);
334 tags[idx][3] = tag_y(3);
335 }
else if (tag_x(0) == 0 && tag_y(0) == 0 && tag_z(0) == 1) {
338 tags[idx][1] = edge_z[tag_x(1)][tag_y(1)];
339 tags[idx][2] = tag_z(2);
340 tags[idx][3] = tag_z(3);
341 }
else if (tag_x(0) == 0 && tag_y(0) == 1 && tag_z(0) == 1) {
344 tags[idx][1] = face_yz[tag_x(1)];
345 tags[idx][2] = tag_y(2) + tag_y(3)*tag_z(2);
346 tags[idx][3] = tag_y(3)*tag_z(3);
347 }
else if (tag_x(0) == 1 && tag_y(0) == 0 && tag_z(0) == 1) {
350 tags[idx][1] = face_xz[tag_y(1)];
351 tags[idx][2] = tag_x(2) + tag_x(3)*tag_z(2);
352 tags[idx][3] = tag_x(3)*tag_z(3);
353 }
else if (tag_x(0) == 1 && tag_y(0) == 1 && tag_z(0) == 0) {
356 tags[idx][1] = face_xy[tag_z(1)];
357 tags[idx][2] = tag_x(2) + tag_x(3)*tag_y(2);
358 tags[idx][3] = tag_x(3)*tag_y(3);
363 tags[idx][2] = tag_x(2) + tag_x(3)*tag_y(2) + tag_x(3)*tag_y(3)*tag_z(2);
364 tags[idx][3] = tag_x(3)*tag_y(3)*tag_z(3);
375 this->setOrdinalTagData(this->tagToOrdinal_,
378 this->basisCardinality_,
386 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
387 dofCoordsHost(
"dofCoordsHost", this->basisCardinality_, this->basisCellTopology_.getDimension());
389 Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>
390 dofCoordsLine(
"dofCoordsLine", cardLine, 1);
392 lineBasis.getDofCoords(dofCoordsLine);
393 auto dofCoordsLineHost = Kokkos::create_mirror_view(dofCoordsLine);
394 Kokkos::deep_copy(dofCoordsLineHost, dofCoordsLine);
396 ordinal_type idx = 0;
397 for (
auto k=0;k<cardLine;++k) {
398 for (ordinal_type j=0;j<cardLine;++j) {
399 for (ordinal_type i=0;i<cardLine;++i,++idx) {
400 dofCoordsHost(idx,0) = dofCoordsLineHost(i,0);
401 dofCoordsHost(idx,1) = dofCoordsLineHost(j,0);
402 dofCoordsHost(idx,2) = dofCoordsLineHost(k,0);
408 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoordsHost);
409 Kokkos::deep_copy(this->dofCoords_, dofCoordsHost);
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Basis_HGRAD_HEX_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
const OrdinalTypeArrayStride1DHost getDofTag(const ordinal_type dofOrd) const
DoF ordinal to DoF tag lookup.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.