16 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
17 #define __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
26 template<EOperator OpType>
27 template<
typename OutputViewType,
28 typename InputViewType,
29 typename WorkViewType,
30 typename VinvViewType>
31 KOKKOS_INLINE_FUNCTION
33 Basis_HGRAD_TRI_Cn_FEM::Serial<OpType>::
34 getValues( OutputViewType output,
35 const InputViewType input,
37 const VinvViewType vinv,
38 const ordinal_type order ) {
40 constexpr ordinal_type spaceDim = 2;
42 card = vinv.extent(0),
43 npts = input.extent(0);
45 typedef typename Kokkos::DynRankView<typename InputViewType::value_type, typename WorkViewType::memory_space> ViewType;
46 auto vcprop = Kokkos::common_view_alloc_prop(input);
47 auto ptr = work.data();
50 case OPERATOR_VALUE: {
51 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
54 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
55 Serial<OpType>::getValues(phis, input, dummyView, order);
57 for (ordinal_type i=0;i<card;++i)
58 for (ordinal_type j=0;j<npts;++j) {
59 output.access(i,j) = 0.0;
60 for (ordinal_type k=0;k<card;++k)
61 output.access(i,j) += vinv(k,i)*phis.access(k,j);
67 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
68 ptr += card*npts*spaceDim*get_dimension_scalar(input);
69 const ViewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
70 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
71 Serial<OpType>::getValues(phis, input, workView, order);
73 for (ordinal_type i=0;i<card;++i)
74 for (ordinal_type j=0;j<npts;++j)
75 for (ordinal_type k=0;k<spaceDim;++k) {
76 output.access(i,j,k) = 0.0;
77 for (ordinal_type l=0;l<card;++l)
78 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
91 const ordinal_type dkcard = getDkCardinality<OpType,spaceDim>();
92 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
95 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
96 Serial<OpType>::getValues(phis, input, dummyView, order);
98 for (ordinal_type i=0;i<card;++i)
99 for (ordinal_type j=0;j<npts;++j)
100 for (ordinal_type k=0;k<dkcard;++k) {
101 output.access(i,j,k) = 0.0;
102 for (ordinal_type l=0;l<card;++l)
103 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
107 case OPERATOR_CURL: {
108 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
109 ptr += card*npts*spaceDim*get_dimension_scalar(input);
110 const ViewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
113 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
114 Serial<OPERATOR_D1>::getValues(phis, input, workView, order);
116 for (ordinal_type i=0;i<card;++i)
117 for (ordinal_type j=0;j<npts;++j) {
118 output.access(i,j,0) = 0.0;
119 for (ordinal_type l=0;l<card;++l)
120 output.access(i,j,0) += vinv(l,i)*phis.access(l,j,1);
121 output.access(i,j,1) = 0.0;
122 for (ordinal_type l=0;l<card;++l)
123 output.access(i,j,1) -= vinv(l,i)*phis.access(l,j,0);
128 INTREPID2_TEST_FOR_ABORT(
true,
129 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented");
134 template<
typename DT, ordinal_type numPtsPerEval,
135 typename outputValueValueType,
class ...outputValueProperties,
136 typename inputPointValueType,
class ...inputPointProperties,
137 typename vinvValueType,
class ...vinvProperties>
139 Basis_HGRAD_TRI_Cn_FEM::
141 const typename DT::execution_space& space,
142 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
143 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
144 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
145 const ordinal_type order,
146 const EOperator operatorType) {
147 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
148 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
149 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
150 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
153 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
154 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
155 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
156 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
158 typedef typename inputPointViewType::value_type inputPointType;
160 const ordinal_type cardinality = outputValues.extent(0);
161 const ordinal_type spaceDim = 2;
163 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
164 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
166 switch (operatorType) {
167 case OPERATOR_VALUE: {
168 workViewType work(Kokkos::view_alloc(space,
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
169 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
170 OPERATOR_VALUE,numPtsPerEval> FunctorType;
171 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
176 workViewType work(Kokkos::view_alloc(space,
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
177 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
178 OPERATOR_D1,numPtsPerEval> FunctorType;
179 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
182 case OPERATOR_CURL: {
183 workViewType work(Kokkos::view_alloc(space,
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
184 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
185 OPERATOR_CURL,numPtsPerEval> FunctorType;
186 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
190 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
191 OPERATOR_D2,numPtsPerEval> FunctorType;
192 workViewType work(Kokkos::view_alloc(space,
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
193 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work, order) );
197 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
198 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented" );
205 template<
typename DT,
typename OT,
typename PT>
208 const EPointType pointType ) {
209 constexpr ordinal_type spaceDim = 2;
211 this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order);
212 this->basisDegree_ = order;
213 this->basisCellTopologyKey_ = shards::Triangle<3>::key;
214 this->basisType_ = BASIS_FEM_LAGRANGIAN;
215 this->basisCoordinates_ = COORDINATES_CARTESIAN;
216 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
218 pointType_ = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
219 const ordinal_type card = this->basisCardinality_;
222 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
223 dofCoords(
"Hgrad::Tri::Cn::dofCoords", card, spaceDim);
226 const shards::CellTopology cellTopo(shards::getCellTopologyData<shards::Triangle<3>>());
227 const ordinal_type offset = 0;
233 this->dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
234 Kokkos::deep_copy(this->dofCoords_, dofCoords);
238 const ordinal_type lwork = card*card;
239 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
240 vmat(
"Hgrad::Tri::Cn::vmat", card, card),
241 work(
"Hgrad::Tri::Cn::work", lwork),
242 ipiv(
"Hgrad::Tri::Cn::ipiv", card);
244 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(
typename Kokkos::HostSpace::execution_space{},
250 ordinal_type info = 0;
251 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
253 lapack.GETRF(card, card,
254 vmat.data(), vmat.stride_1(),
255 (ordinal_type*)ipiv.data(),
258 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
260 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
263 vmat.data(), vmat.stride_1(),
264 (ordinal_type*)ipiv.data(),
268 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
270 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
273 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
274 vinv(
"Hgrad::Line::Cn::vinv", card, card);
276 for (ordinal_type i=0;i<card;++i)
277 for (ordinal_type j=0;j<card;++j)
278 vinv(i,j) = vmat(j,i);
280 this->vinv_ = Kokkos::create_mirror_view(
typename DT::memory_space(), vinv);
281 Kokkos::deep_copy(this->vinv_ , vinv);
286 constexpr ordinal_type tagSize = 4;
287 const ordinal_type posScDim = 0;
288 const ordinal_type posScOrd = 1;
289 const ordinal_type posDfOrd = 2;
293 INTREPID2_TEST_FOR_EXCEPTION( order >
Parameters::MaxOrder, std::invalid_argument,
"polynomial order exceeds the max supported by this class");
294 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
295 ordinal_type tags[maxCard][tagSize];
298 numEdgeDof = Intrepid2::getPnCardinality<1>(order-2),
299 numElemDof = (order > 2 ? Intrepid2::getPnCardinality<2>(order-3) : 0);
304 ordinal_type edgeId[3] = {}, elemId = 0;
305 for (ordinal_type i=0;i<card;++i) {
308 const auto x = dofCoords(i,0);
309 const auto y = dofCoords(i,1);
315 if ((1.0 - xi0) < eps) {
321 else if ((1.0 - xi1) < eps) {
327 else if ((1.0 - xi2) < eps) {
333 else if (xi2 < eps) {
336 tags[i][2] = edgeId[0]++;
337 tags[i][3] = numEdgeDof;
339 else if (xi0 < eps) {
342 tags[i][2] = edgeId[1]++;
343 tags[i][3] = numEdgeDof;
345 else if (xi1 < eps) {
348 tags[i][2] = edgeId[2]++;
349 tags[i][3] = numEdgeDof;
354 tags[i][2] = elemId++;
355 tags[i][3] = numElemDof;
363 this->setOrdinalTagData(this->tagToOrdinal_,
366 this->basisCardinality_,
374 template<
typename DT,
typename OT,
typename PT>
377 ordinal_type& perTeamSpaceSize,
378 ordinal_type& perThreadSpaceSize,
380 const EOperator operatorType)
const {
381 perTeamSpaceSize = 0;
382 perThreadSpaceSize = getWorkSizePerPoint(operatorType)*get_dimension_scalar(inputPoints)*
sizeof(
typename BasisBase::scalarType);
385 template<
typename DT,
typename OT,
typename PT>
386 KOKKOS_INLINE_FUNCTION
389 OutputViewType outputValues,
390 const PointViewType inputPoints,
391 const EOperator operatorType,
392 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
393 const typename DT::execution_space::scratch_memory_space & scratchStorage,
394 const ordinal_type subcellDim,
395 const ordinal_type subcellOrdinal)
const {
397 INTREPID2_TEST_FOR_ABORT( !((subcellDim == -1) && (subcellOrdinal == -1)),
398 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
400 const int numPoints = inputPoints.extent(0);
401 using ScalarType =
typename ScalarTraits<typename PointViewType::value_type>::scalar_type;
402 using WorkViewType = Kokkos::DynRankView< ScalarType,typename DT::execution_space::scratch_memory_space,Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
403 constexpr ordinal_type spaceDim = 2;
404 auto sizePerPoint = (operatorType==OPERATOR_VALUE) ?
405 this->vinv_.extent(0)*get_dimension_scalar(inputPoints) :
406 (2*spaceDim+1)*this->vinv_.extent(0)*get_dimension_scalar(inputPoints);
407 WorkViewType workView(scratchStorage, sizePerPoint*team_member.team_size());
408 using range_type = Kokkos::pair<ordinal_type,ordinal_type>;
409 switch(operatorType) {
411 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
412 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type (pt,pt+1), Kokkos::ALL() );
413 const auto input = Kokkos::subview( inputPoints, range_type(pt, pt+1), Kokkos::ALL() );
414 WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
419 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
420 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type(pt,pt+1), Kokkos::ALL() );
421 const auto input = Kokkos::subview( inputPoints, range_type(pt,pt+1), Kokkos::ALL() );
422 WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
423 Impl::Basis_HGRAD_TRI_Cn_FEM::Serial<OPERATOR_GRAD>::getValues( output, input, work, this->vinv_, this->basisDegree_);
427 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
428 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type(pt,pt+1), Kokkos::ALL() );
429 const auto input = Kokkos::subview( inputPoints, range_type(pt,pt+1), Kokkos::ALL() );
430 WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
431 Impl::Basis_HGRAD_TRI_Cn_FEM::Serial<OPERATOR_CURL>::getValues( output, input, work, this->vinv_, this->basisDegree_);
435 INTREPID2_TEST_FOR_ABORT(
true,
436 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): getValues not implemented for this operator");
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Return the size of the scratch space, in bytes, needed for using the team-level implementation of get...
See Intrepid2::Basis_HGRAD_TRI_Cn_FEM work is a rank 1 view having the same value_type of inputPoints...
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.