15 #ifndef __INTREPID2_HVOL_LINE_CN_FEM_DEF_HPP__
16 #define __INTREPID2_HVOL_LINE_CN_FEM_DEF_HPP__
23 template<EOperator opType>
24 template<
typename OutputViewType,
25 typename InputViewType,
26 typename WorkViewType,
27 typename VinvViewType>
28 KOKKOS_INLINE_FUNCTION
30 Basis_HVOL_LINE_Cn_FEM::Serial<opType>::
31 getValues( OutputViewType output,
32 const InputViewType input,
34 const VinvViewType vinv,
35 const ordinal_type operatorDn ) {
36 ordinal_type opDn = operatorDn;
38 const ordinal_type card = vinv.extent(0);
39 const ordinal_type npts = input.extent(0);
41 const ordinal_type order = card - 1;
42 const double alpha = 0.0, beta = 0.0;
44 typedef typename Kokkos::DynRankView<typename InputViewType::value_type, typename WorkViewType::memory_space> ViewType;
45 auto vcprop = Kokkos::common_view_alloc_prop(input);
48 case OPERATOR_VALUE: {
49 ViewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts);
51 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
52 Serial<opType>::getValues(phis, input, order, alpha, beta);
54 for (ordinal_type i=0;i<card;++i)
55 for (ordinal_type j=0;j<npts;++j) {
56 output.access(i,j) = 0.0;
57 for (ordinal_type k=0;k<card;++k)
58 output.access(i,j) += vinv(k,i)*phis.access(k,j);
73 opDn = getOperatorOrder(opType);
76 const ordinal_type dkcard = 1;
77 ViewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts, dkcard);
78 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
79 Serial<opType>::getValues(phis, input, order, alpha, beta, opDn);
81 for (ordinal_type i=0;i<card;++i)
82 for (ordinal_type j=0;j<npts;++j)
83 for (ordinal_type k=0;k<dkcard;++k) {
84 output.access(i,j,k) = 0.0;
85 for (ordinal_type l=0;l<card;++l)
86 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
91 INTREPID2_TEST_FOR_ABORT(
true,
92 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::Serial::getValues) operator is not supported." );
98 template<
typename DT, ordinal_type numPtsPerEval,
99 typename outputValueValueType,
class ...outputValueProperties,
100 typename inputPointValueType,
class ...inputPointProperties,
101 typename vinvValueType,
class ...vinvProperties>
103 Basis_HVOL_LINE_Cn_FEM::
104 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
105 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
106 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
107 const EOperator operatorType ) {
108 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
109 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
110 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
111 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
114 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
115 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
116 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
117 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
119 typedef typename inputPointViewType::value_type inputPointType;
121 const ordinal_type cardinality = outputValues.extent(0);
123 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
124 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
125 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_LINE_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
127 switch (operatorType) {
128 case OPERATOR_VALUE: {
129 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
130 OPERATOR_VALUE,numPtsPerEval> FunctorType;
131 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
145 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
146 OPERATOR_Dn,numPtsPerEval> FunctorType;
147 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
148 getOperatorOrder(operatorType)) );
152 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
153 ">>> ERROR (Basis_HVOL_LINE_Cn_FEM): Operator type not implemented" );
162 template<
typename DT,
typename OT,
typename PT>
165 const EPointType pointType ) {
166 this->pointType_ = pointType;
167 this->basisCardinality_ = order+1;
168 this->basisDegree_ = order;
169 this->basisCellTopologyKey_ = shards::Line<2>::key;
170 this->basisType_ = BASIS_FEM_LAGRANGIAN;
171 this->basisCoordinates_ = COORDINATES_CARTESIAN;
172 this->functionSpace_ = FUNCTION_SPACE_HVOL;
174 const ordinal_type card = this->basisCardinality_;
177 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
178 dofCoords(
"HVOL::Line::Cn::dofCoords", card, 1);
181 auto pointT = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
184 case POINTTYPE_EQUISPACED:
185 case POINTTYPE_WARPBLEND: {
188 const shards::CellTopology cellTopo(shards::getCellTopologyData<shards::Line<2> >());
189 const ordinal_type offset = 1;
192 order+1+offset, offset,
198 case POINTTYPE_GAUSS: {
205 INTREPID2_TEST_FOR_EXCEPTION( !isValidPointType(pointT),
206 std::invalid_argument ,
207 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM) invalid pointType." );
211 this->dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
212 Kokkos::deep_copy(this->dofCoords_, dofCoords);
216 const ordinal_type lwork = card*card;
217 Kokkos::DynRankView<typename ScalarViewType::value_type,Kokkos::LayoutLeft,Kokkos::HostSpace>
218 vmat(
"HVOL::Line::Cn::vmat", card, card),
219 work(
"HVOL::Line::Cn::work", lwork),
220 ipiv(
"HVOL::Line::Cn::ipiv", card);
222 const double alpha = 0.0, beta = 0.0;
223 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
224 getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>
225 (
typename Kokkos::HostSpace::execution_space{}, vmat, dofCoords, order, alpha, beta, OPERATOR_VALUE);
227 ordinal_type info = 0;
228 Teuchos::LAPACK<ordinal_type,typename ScalarViewType::value_type> lapack;
230 lapack.GETRF(card, card,
231 vmat.data(), vmat.stride_1(),
232 (ordinal_type*)ipiv.data(),
235 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
237 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM) lapack.GETRF returns nonzero info." );
240 vmat.data(), vmat.stride_1(),
241 (ordinal_type*)ipiv.data(),
245 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
247 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM) lapack.GETRI returns nonzero info." );
250 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
251 vinv(
"HVOL::Line::Cn::vinv", card, card);
253 for (ordinal_type i=0;i<card;++i)
254 for (ordinal_type j=0;j<card;++j)
255 vinv(i,j) = vmat(j,i);
257 this->vinv_ = Kokkos::create_mirror_view(
typename DT::memory_space(), vinv);
258 Kokkos::deep_copy(this->vinv_ , vinv);
263 const ordinal_type tagSize = 4;
264 const ordinal_type posScDim = 0;
265 const ordinal_type posScOrd = 1;
266 const ordinal_type posDfOrd = 2;
271 for (ordinal_type i=0;i<card;++i) {
282 this->setOrdinalTagData(this->tagToOrdinal_,
285 this->basisCardinality_,
293 template<
typename DT,
typename OT,
typename PT>
296 ordinal_type& perTeamSpaceSize,
297 ordinal_type& perThreadSpaceSize,
299 const EOperator operatorType)
const {
300 perTeamSpaceSize = 0;
301 perThreadSpaceSize = this->vinv_.extent(0)*get_dimension_scalar(inputPoints)*
sizeof(
typename BasisBase::scalarType);
304 template<
typename DT,
typename OT,
typename PT>
305 KOKKOS_INLINE_FUNCTION
308 OutputViewType outputValues,
309 const PointViewType inputPoints,
310 const EOperator operatorType,
311 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
312 const typename DT::execution_space::scratch_memory_space & scratchStorage,
313 const ordinal_type subcellDim,
314 const ordinal_type subcellOrdinal)
const {
316 INTREPID2_TEST_FOR_ABORT( !((subcellDim == -1) && (subcellOrdinal == -1)),
317 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
319 const int numPoints = inputPoints.extent(0);
320 using ScalarType =
typename ScalarTraits<typename PointViewType::value_type>::scalar_type;
321 using WorkViewType = Kokkos::DynRankView< ScalarType,typename DT::execution_space::scratch_memory_space,Kokkos::MemoryTraits<Kokkos::Unmanaged> >;
322 ordinal_type sizePerPoint = this->vinv_.extent(0)*get_dimension_scalar(inputPoints);
323 WorkViewType workView(scratchStorage, sizePerPoint*team_member.team_size());
324 using range_type = Kokkos::pair<ordinal_type,ordinal_type>;
326 switch(operatorType) {
328 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
329 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), range_type (pt,pt+1), Kokkos::ALL() );
330 const auto input = Kokkos::subview( inputPoints, range_type(pt, pt+1), Kokkos::ALL() );
331 WorkViewType work(workView.data() + sizePerPoint*team_member.team_rank(), sizePerPoint);
336 INTREPID2_TEST_FOR_ABORT(
true,
337 ">>> ERROR (Basis_HVOL_LINE_Cn_FEM): getValues not implemented for this operator");
ScalarTraits< pointValueType >::scalar_type scalarType
Scalar type for point values.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
virtual void getScratchSpaceSize(ordinal_type &perTeamSpaceSize, ordinal_type &perThreadSpaceSize, const PointViewType inputPointsconst, 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...
virtual void getValues(OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
See Intrepid2::Basis_HVOL_LINE_Cn_FEM.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.