15 #ifndef __INTREPID2_HVOL_TRI_CN_FEM_DEF_HPP__
16 #define __INTREPID2_HVOL_TRI_CN_FEM_DEF_HPP__
25 template<EOperator OpType>
26 template<
typename OutputViewType,
27 typename InputViewType,
28 typename WorkViewType,
29 typename VinvViewType>
30 KOKKOS_INLINE_FUNCTION
32 Basis_HVOL_TRI_Cn_FEM::Serial<OpType>::
33 getValues( OutputViewType output,
34 const InputViewType input,
36 const VinvViewType vinv ) {
38 constexpr ordinal_type spaceDim = 2;
40 card = vinv.extent(0),
41 npts = input.extent(0);
44 ordinal_type order = 0;
46 if (card == Intrepid2::getPnCardinality<spaceDim>(p) ) {
52 typedef typename Kokkos::DynRankView<typename InputViewType::value_type, typename WorkViewType::memory_space> ViewType;
53 auto vcprop = Kokkos::common_view_alloc_prop(input);
54 auto ptr = work.data();
57 case OPERATOR_VALUE: {
58 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
61 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
62 Serial<OpType>::getValues(phis, input, dummyView, order);
64 for (ordinal_type i=0;i<card;++i)
65 for (ordinal_type j=0;j<npts;++j) {
66 output.access(i,j) = 0.0;
67 for (ordinal_type k=0;k<card;++k)
68 output.access(i,j) += vinv(k,i)*phis.access(k,j);
74 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
75 ptr += card*npts*spaceDim*get_dimension_scalar(input);
76 const ViewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
77 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
78 Serial<OpType>::getValues(phis, input, workView, order);
80 for (ordinal_type i=0;i<card;++i)
81 for (ordinal_type j=0;j<npts;++j)
82 for (ordinal_type k=0;k<spaceDim;++k) {
83 output.access(i,j,k) = 0.0;
84 for (ordinal_type l=0;l<card;++l)
85 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
98 const ordinal_type dkcard = getDkCardinality<OpType,spaceDim>();
99 const ViewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
102 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
103 Serial<OpType>::getValues(phis, input, dummyView, order);
105 for (ordinal_type i=0;i<card;++i)
106 for (ordinal_type j=0;j<npts;++j)
107 for (ordinal_type k=0;k<dkcard;++k) {
108 output.access(i,j,k) = 0.0;
109 for (ordinal_type l=0;l<card;++l)
110 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
115 INTREPID2_TEST_FOR_ABORT(
true,
116 ">>> ERROR (Basis_HVOL_TRI_Cn_FEM): Operator type not implemented");
121 template<
typename DT, ordinal_type numPtsPerEval,
122 typename outputValueValueType,
class ...outputValueProperties,
123 typename inputPointValueType,
class ...inputPointProperties,
124 typename vinvValueType,
class ...vinvProperties>
126 Basis_HVOL_TRI_Cn_FEM::
127 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
130 const EOperator operatorType) {
131 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
132 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
133 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
134 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
137 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
138 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
139 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
140 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
142 typedef typename inputPointViewType::value_type inputPointType;
144 const ordinal_type cardinality = outputValues.extent(0);
145 const ordinal_type spaceDim = 2;
147 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
148 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
150 switch (operatorType) {
151 case OPERATOR_VALUE: {
152 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
153 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
154 OPERATOR_VALUE,numPtsPerEval> FunctorType;
155 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
160 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
161 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
162 OPERATOR_D1,numPtsPerEval> FunctorType;
163 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
167 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
168 OPERATOR_D2,numPtsPerEval> FunctorType;
169 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
170 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
181 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
182 ">>> ERROR (Basis_HVOL_TRI_Cn_FEM): Operator type not implemented" );
189 template<
typename DT,
typename OT,
typename PT>
192 const EPointType pointType ) {
194 constexpr ordinal_type spaceDim = 2;
196 this->pointType_ = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
197 this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order);
198 this->basisDegree_ = order;
199 this->basisCellTopologyKey_ = shards::Triangle<3>::key;
200 this->basisType_ = BASIS_FEM_LAGRANGIAN;
201 this->basisCoordinates_ = COORDINATES_CARTESIAN;
202 this->functionSpace_ = FUNCTION_SPACE_HVOL;
204 const ordinal_type card = this->basisCardinality_;
207 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
208 dofCoords(
"HVOL::Tri::Cn::dofCoords", card, spaceDim);
211 const ordinal_type offset = 1;
212 const shards::CellTopology cellTopo(shards::getCellTopologyData<shards::Triangle<3> >());
215 order+spaceDim+offset, offset,
218 this->dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
219 Kokkos::deep_copy(this->dofCoords_, dofCoords);
223 const ordinal_type lwork = card*card;
224 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
225 vmat(
"HVOL::Tri::Cn::vmat", card, card),
226 work(
"HVOL::Tri::Cn::work", lwork),
227 ipiv(
"HVOL::Tri::Cn::ipiv", card);
229 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(
typename Kokkos::HostSpace::execution_space{},
235 ordinal_type info = 0;
236 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
238 lapack.GETRF(card, card,
239 vmat.data(), vmat.stride_1(),
240 (ordinal_type*)ipiv.data(),
243 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
245 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
248 vmat.data(), vmat.stride_1(),
249 (ordinal_type*)ipiv.data(),
253 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
255 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
258 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
259 vinv(
"HVOL::Line::Cn::vinv", card, card);
261 for (ordinal_type i=0;i<card;++i)
262 for (ordinal_type j=0;j<card;++j)
263 vinv(i,j) = vmat(j,i);
265 this->vinv_ = Kokkos::create_mirror_view(
typename DT::memory_space(), vinv);
266 Kokkos::deep_copy(this->vinv_ , vinv);
271 constexpr ordinal_type tagSize = 4;
272 const ordinal_type posScDim = 0;
273 const ordinal_type posScOrd = 1;
274 const ordinal_type posDfOrd = 2;
276 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
277 ordinal_type tags[maxCard][tagSize];
280 numElemDof = this->basisCardinality_;
282 ordinal_type elemId = 0;
283 for (ordinal_type i=0;i<this->basisCardinality_;++i) {
285 tags[i][0] = spaceDim;
287 tags[i][2] = elemId++;
288 tags[i][3] = numElemDof;
295 this->setOrdinalTagData(this->tagToOrdinal_,
298 this->basisCardinality_,
306 template<
typename DT,
typename OT,
typename PT>
309 ordinal_type& perTeamSpaceSize,
310 ordinal_type& perThreadSpaceSize,
312 const EOperator operatorType)
const {
313 perTeamSpaceSize = 0;
314 perThreadSpaceSize = this->vinv_.extent(0)*get_dimension_scalar(inputPoints)*
sizeof(
typename BasisBase::scalarType);
317 template<
typename DT,
typename OT,
typename PT>
318 KOKKOS_INLINE_FUNCTION
321 OutputViewType outputValues,
322 const PointViewType inputPoints,
323 const EOperator operatorType,
324 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
325 const typename DT::execution_space::scratch_memory_space & scratchStorage,
326 const ordinal_type subcellDim,
327 const ordinal_type subcellOrdinal)
const {
329 INTREPID2_TEST_FOR_ABORT( !((subcellDim == -1) && (subcellOrdinal == -1)),
330 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
332 const int numPoints = inputPoints.extent(0);
333 using ScalarType =
typename ScalarTraits<typename PointViewType::value_type>::scalar_type;
334 using WorkViewType = Kokkos::DynRankView< ScalarType,typename DT::execution_space::scratch_memory_space,Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
335 auto sizePerPoint = this->vinv_.extent(0)*get_dimension_scalar(inputPoints);
336 WorkViewType workView(scratchStorage, sizePerPoint*team_member.team_size());
337 using range_type = Kokkos::pair<ordinal_type,ordinal_type>;
338 switch(operatorType) {
340 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
341 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type (pt,pt+1), Kokkos::ALL() );
342 const auto input = Kokkos::subview( inputPoints, range_type(pt, pt+1), Kokkos::ALL() );
343 WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
348 INTREPID2_TEST_FOR_ABORT(
true,
349 ">>> ERROR (Basis_HVOL_TRI_Cn_FEM): getValues not implemented for this operator");
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
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...
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
See Intrepid2::Basis_HVOL_TRI_Cn_FEM.