48 #ifndef __INTREPID2_HVOL_HEX_CN_FEMDEF_HPP__
49 #define __INTREPID2_HVOL_HEX_CN_FEMDEF_HPP__
56 template<EOperator opType>
57 template<
typename OutputViewType,
58 typename inputViewType,
59 typename workViewType,
60 typename vinvViewType>
61 KOKKOS_INLINE_FUNCTION
63 Basis_HVOL_HEX_Cn_FEM::Serial<opType>::
64 getValues( OutputViewType output,
65 const inputViewType input,
67 const vinvViewType vinv,
68 const ordinal_type operatorDn ) {
69 ordinal_type opDn = operatorDn;
71 const ordinal_type cardLine = vinv.extent(0);
72 const ordinal_type npts = input.extent(0);
74 typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
75 const auto input_x = Kokkos::subview(input, Kokkos::ALL(), range_type(0,1));
76 const auto input_y = Kokkos::subview(input, Kokkos::ALL(), range_type(1,2));
77 const auto input_z = Kokkos::subview(input, Kokkos::ALL(), range_type(2,3));
79 const ordinal_type dim_s = get_dimension_scalar(work);
80 auto ptr0 = work.data();
81 auto ptr1 = work.data()+cardLine*npts*dim_s;
82 auto ptr2 = work.data()+2*cardLine*npts*dim_s;
83 auto ptr3 = work.data()+3*cardLine*npts*dim_s;
85 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
86 auto vcprop = Kokkos::common_view_alloc_prop(work);
89 case OPERATOR_VALUE: {
90 viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
91 viewType output_x(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
92 viewType output_y(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
93 viewType output_z(Kokkos::view_wrap(ptr3, vcprop), cardLine, npts);
95 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
96 getValues(output_x, input_x, work_line, vinv);
98 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
99 getValues(output_y, input_y, work_line, vinv);
101 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
102 getValues(output_z, input_z, work_line, vinv);
105 ordinal_type idx = 0;
106 for (ordinal_type k=0;k<cardLine;++k)
107 for (ordinal_type j=0;j<cardLine;++j)
108 for (ordinal_type i=0;i<cardLine;++i,++idx)
109 for (ordinal_type l=0;l<npts;++l)
110 output.access(idx,l) = output_x.access(i,l)*output_y.access(j,l)*output_z.access(k,l);
124 opDn = getOperatorOrder(opType);
126 const ordinal_type dkcard = opDn + 1;
129 for (ordinal_type l1=0;l1<dkcard;++l1)
130 for (ordinal_type l0=0;l0<(l1+1);++l0) {
131 const ordinal_type mult_x = (opDn - l1);
132 const ordinal_type mult_y = l1 - l0;
133 const ordinal_type mult_z = l0;
141 viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
142 decltype(work_line) output_x, output_y, output_z;
145 output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts, 1);
146 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
147 getValues(output_x, input_x, work_line, vinv, mult_x);
149 output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
150 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
151 getValues(output_x, input_x, work_line, vinv);
155 output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts, 1);
156 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
157 getValues(output_y, input_y, work_line, vinv, mult_y);
159 output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
160 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
161 getValues(output_y, input_y, work_line, vinv);
165 output_z = viewType(Kokkos::view_wrap(ptr3, vcprop), cardLine, npts, 1);
166 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
167 getValues(output_z, input_z, work_line, vinv, mult_z);
169 output_z = viewType(Kokkos::view_wrap(ptr3, vcprop), cardLine, npts);
170 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
171 getValues(output_z, input_z, work_line, vinv);
175 ordinal_type idx = 0;
176 for (ordinal_type k=0;k<cardLine;++k)
177 for (ordinal_type j=0;j<cardLine;++j)
178 for (ordinal_type i=0;i<cardLine;++i,++idx)
179 for (ordinal_type l=0;l<npts;++l)
180 output.access(idx,l,d) = output_x.access(i,l,0)*output_y.access(j,l,0)*output_z.access(k,l,0);
187 INTREPID2_TEST_FOR_ABORT(
true ,
188 ">>> ERROR (Basis_HVOL_HEX_Cn_FEM): Operator type not implemented");
194 template<
typename SpT, ordinal_type numPtsPerEval,
195 typename outputValueValueType,
class ...outputValueProperties,
196 typename inputPointValueType,
class ...inputPointProperties,
197 typename vinvValueType,
class ...vinvProperties>
199 Basis_HVOL_HEX_Cn_FEM::
200 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
201 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
202 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
203 const EOperator operatorType ) {
204 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
205 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
206 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
207 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
210 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
211 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
212 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
213 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
215 typedef typename inputPointViewType::value_type inputPointType;
217 const ordinal_type cardinality = outputValues.extent(0);
218 const ordinal_type cardLine = std::cbrt(cardinality);
219 const ordinal_type workSize = 4*cardLine;
221 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
222 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
223 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_HEX_Cn_FEM::getValues::work", vcprop), workSize, inputPoints.extent(0));
225 switch (operatorType) {
226 case OPERATOR_VALUE: {
227 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
228 OPERATOR_VALUE,numPtsPerEval> FunctorType;
229 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
243 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
244 OPERATOR_Dn,numPtsPerEval> FunctorType;
245 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
246 getOperatorOrder(operatorType)) );
250 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
251 ">>> ERROR (Basis_HVOL_HEX_Cn_FEM): Operator type not implemented" );
259 template<
typename SpT,
typename OT,
typename PT>
262 const EPointType pointType ) {
268 this->vinv_ = Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>(
"HVOL::HEX::Cn::vinv", cardLine, cardLine);
269 lineBasis.getVandermondeInverse(this->vinv_);
271 this->basisCardinality_ = cardLine*cardLine*cardLine;
272 this->basisDegree_ = order;
273 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Hexahedron<8> >() );
274 this->basisType_ = BASIS_FEM_FIAT;
275 this->basisCoordinates_ = COORDINATES_CARTESIAN;
276 this->functionSpace_ = FUNCTION_SPACE_HVOL;
281 const ordinal_type tagSize = 4;
282 const ordinal_type posScDim = 0;
283 const ordinal_type posScOrd = 1;
284 const ordinal_type posDfOrd = 2;
288 ordinal_type tags[maxCardLine*maxCardLine*maxCardLine][4];
291 ordinal_type idx = 0;
292 for (
auto k=0;k<cardLine;++k) {
293 const auto tag_z = lineBasis.
getDofTag(k);
294 for (ordinal_type j=0;j<cardLine;++j) {
295 const auto tag_y = lineBasis.
getDofTag(j);
296 for (ordinal_type i=0;i<cardLine;++i,++idx) {
297 const auto tag_x = lineBasis.
getDofTag(i);
302 tags[idx][2] = tag_x(2) + tag_x(3)*tag_y(2) + tag_x(3)*tag_y(3)*tag_z(2);
303 tags[idx][3] = tag_x(3)*tag_y(3)*tag_z(3);
313 this->setOrdinalTagData(this->tagToOrdinal_,
316 this->basisCardinality_,
324 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
325 dofCoordsHost(
"dofCoordsHost", this->basisCardinality_, this->basisCellTopology_.getDimension());
327 Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>
328 dofCoordsLine(
"dofCoordsLine", cardLine, 1);
330 lineBasis.getDofCoords(dofCoordsLine);
331 auto dofCoordsLineHost = Kokkos::create_mirror_view(Kokkos::HostSpace(), dofCoordsLine);
332 Kokkos::deep_copy(dofCoordsLineHost, dofCoordsLine);
334 ordinal_type idx = 0;
335 for (
auto k=0;k<cardLine;++k) {
336 for (ordinal_type j=0;j<cardLine;++j) {
337 for (ordinal_type i=0;i<cardLine;++i,++idx) {
338 dofCoordsHost(idx,0) = dofCoordsLineHost(i,0);
339 dofCoordsHost(idx,1) = dofCoordsLineHost(j,0);
340 dofCoordsHost(idx,2) = dofCoordsLineHost(k,0);
347 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoordsHost);
348 Kokkos::deep_copy(this->dofCoords_, dofCoordsHost);
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Basis_HVOL_HEX_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
const OrdinalTypeArrayStride1DHost getDofTag(const ordinal_type dofOrd) const
DoF ordinal to DoF tag lookup.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.