48 #ifndef __INTREPID2_HVOL_QUAD_CN_FEM_DEF_HPP__
49 #define __INTREPID2_HVOL_QUAD_CN_FEM_DEF_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_QUAD_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));
78 const int dim_s = get_dimension_scalar(work);
79 auto ptr0 = work.data();
80 auto ptr1 = work.data()+cardLine*npts*dim_s;
81 auto ptr2 = work.data()+2*cardLine*npts*dim_s;
83 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
84 auto vcprop = Kokkos::common_view_alloc_prop(work);
87 case OPERATOR_VALUE: {
88 viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
89 viewType output_x(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
90 viewType output_y(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
92 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
93 getValues(output_x, input_x, work_line, vinv);
95 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
96 getValues(output_y, input_y, work_line, vinv);
100 for (ordinal_type j=0;j<cardLine;++j)
101 for (ordinal_type i=0;i<cardLine;++i,++idx)
102 for (ordinal_type k=0;k<npts;++k)
103 output.access(idx,k) = output_x.access(i,k)*output_y.access(j,k);
117 opDn = getOperatorOrder(opType);
119 const auto dkcard = opDn + 1;
120 for (
auto l=0;l<dkcard;++l) {
121 viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
123 viewType output_x, output_y;
125 const auto mult_x = opDn - l;
126 const auto mult_y = l;
129 output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts, 1);
130 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
131 getValues(output_x, input_x, work_line, vinv, mult_x);
133 output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
134 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
135 getValues(output_x, input_x, work_line, vinv);
139 output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts, 1);
140 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
141 getValues(output_y, input_y, work_line, vinv, mult_y);
143 output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
144 Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
145 getValues(output_y, input_y, work_line, vinv);
149 ordinal_type idx = 0;
150 for (ordinal_type j=0;j<cardLine;++j)
151 for (ordinal_type i=0;i<cardLine;++i,++idx)
152 for (ordinal_type k=0;k<npts;++k)
153 output.access(idx,k,l) = output_x.access(i,k,0)*output_y.access(j,k,0);
158 INTREPID2_TEST_FOR_ABORT(
true,
159 ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Serial::getValues) operator is not supported" );
164 template<
typename SpT, ordinal_type numPtsPerEval,
165 typename outputValueValueType,
class ...outputValueProperties,
166 typename inputPointValueType,
class ...inputPointProperties,
167 typename vinvValueType,
class ...vinvProperties>
169 Basis_HVOL_QUAD_Cn_FEM::
170 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
171 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
172 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
173 const EOperator operatorType ) {
174 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
175 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
176 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
177 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
180 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
181 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
182 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
183 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
185 typedef typename inputPointViewType::value_type inputPointType;
187 const ordinal_type cardinality = outputValues.extent(0);
188 const ordinal_type cardLine = std::sqrt(cardinality);
189 const ordinal_type workSize = 3*cardLine;
191 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
192 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
193 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_QUAD_Cn_FEM::getValues::work", vcprop), workSize, inputPoints.extent(0));
195 switch (operatorType) {
196 case OPERATOR_VALUE: {
197 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
198 OPERATOR_VALUE,numPtsPerEval> FunctorType;
199 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
213 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
214 OPERATOR_Dn,numPtsPerEval> FunctorType;
215 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
216 getOperatorOrder(operatorType)) );
220 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
221 ">>> ERROR (Basis_HVOL_QUAD_Cn_FEM): Operator type not implemented" );
229 template<
typename SpT,
typename OT,
typename PT>
232 const EPointType pointType ) {
241 this->vinv_ = Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>(
"HVOL::Quad::Cn::vinv", cardLine, cardLine);
242 lineBasis.getVandermondeInverse(this->vinv_);
244 this->basisCardinality_ = cardLine*cardLine;
245 this->basisDegree_ = order;
246 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
247 this->basisType_ = BASIS_FEM_FIAT;
248 this->basisCoordinates_ = COORDINATES_CARTESIAN;
249 this->functionSpace_ = FUNCTION_SPACE_HVOL;
254 const ordinal_type tagSize = 4;
255 const ordinal_type posScDim = 0;
256 const ordinal_type posScOrd = 1;
257 const ordinal_type posDfOrd = 2;
261 ordinal_type tags[maxCardLine*maxCardLine][4];
264 ordinal_type idx = 0;
265 for (ordinal_type j=0;j<cardLine;++j) {
266 const auto tag_y = lineBasis.
getDofTag(j);
267 for (ordinal_type i=0;i<cardLine;++i,++idx) {
268 const auto tag_x = lineBasis.
getDofTag(i);
273 tags[idx][2] = tag_x(2) + tag_x(3)*tag_y(2);
274 tags[idx][3] = tag_x(3)*tag_y(3);
283 this->setOrdinalTagData(this->tagToOrdinal_,
286 this->basisCardinality_,
294 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
295 dofCoordsHost(
"dofCoordsHost", this->basisCardinality_, this->basisCellTopology_.getDimension());
297 Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>
298 dofCoordsLine(
"dofCoordsLine", cardLine, 1);
300 lineBasis.getDofCoords(dofCoordsLine);
301 auto dofCoordsLineHost = Kokkos::create_mirror_view(dofCoordsLine);
302 Kokkos::deep_copy(dofCoordsLineHost, dofCoordsLine);
304 ordinal_type idx = 0;
305 for (ordinal_type j=0;j<cardLine;++j) {
306 for (ordinal_type i=0;i<cardLine;++i,++idx) {
307 dofCoordsHost(idx,0) = dofCoordsLineHost(i,0);
308 dofCoordsHost(idx,1) = dofCoordsLineHost(j,0);
313 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoordsHost);
314 Kokkos::deep_copy(this->dofCoords_, dofCoordsHost);
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Basis_HVOL_QUAD_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 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.