48 #ifndef __INTREPID2_HVOL_LINE_CN_FEM_DEF_HPP__
49 #define __INTREPID2_HVOL_LINE_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_LINE_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 card = vinv.extent(0);
72 const ordinal_type npts = input.extent(0);
74 const ordinal_type order = card - 1;
75 const double alpha = 0.0, beta = 0.0;
77 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
78 auto vcprop = Kokkos::common_view_alloc_prop(work);
81 case OPERATOR_VALUE: {
82 viewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts);
84 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
85 Serial<opType>::getValues(phis, input, order, alpha, beta);
87 for (ordinal_type i=0;i<card;++i)
88 for (ordinal_type j=0;j<npts;++j) {
89 output.access(i,j) = 0.0;
90 for (ordinal_type k=0;k<card;++k)
91 output.access(i,j) += vinv(k,i)*phis.access(k,j);
106 opDn = getOperatorOrder(opType);
109 const ordinal_type dkcard = 1;
110 viewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts, dkcard);
111 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
112 Serial<opType>::getValues(phis, input, order, alpha, beta, opDn);
114 for (ordinal_type i=0;i<card;++i)
115 for (ordinal_type j=0;j<npts;++j)
116 for (ordinal_type k=0;k<dkcard;++k) {
117 output.access(i,j,k) = 0.0;
118 for (ordinal_type l=0;l<card;++l)
119 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
124 INTREPID2_TEST_FOR_ABORT(
true,
125 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM::Serial::getValues) operator is not supported." );
131 template<
typename SpT, ordinal_type numPtsPerEval,
132 typename outputValueValueType,
class ...outputValueProperties,
133 typename inputPointValueType,
class ...inputPointProperties,
134 typename vinvValueType,
class ...vinvProperties>
136 Basis_HVOL_LINE_Cn_FEM::
137 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
138 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
139 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
140 const EOperator operatorType ) {
141 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
142 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
143 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
144 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
147 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
148 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
149 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
150 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
152 typedef typename inputPointViewType::value_type inputPointType;
154 const ordinal_type cardinality = outputValues.extent(0);
156 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
157 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
158 workViewType work(Kokkos::view_alloc(
"Basis_HVOL_LINE_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
160 switch (operatorType) {
161 case OPERATOR_VALUE: {
162 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
163 OPERATOR_VALUE,numPtsPerEval> FunctorType;
164 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
178 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
179 OPERATOR_Dn,numPtsPerEval> FunctorType;
180 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
181 getOperatorOrder(operatorType)) );
185 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
186 ">>> ERROR (Basis_HVOL_LINE_Cn_FEM): Operator type not implemented" );
195 template<
typename SpT,
typename OT,
typename PT>
198 const EPointType pointType ) {
199 this->basisCardinality_ = order+1;
200 this->basisDegree_ = order;
201 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Line<2> >() );
202 this->basisType_ = BASIS_FEM_FIAT;
203 this->basisCoordinates_ = COORDINATES_CARTESIAN;
204 this->functionSpace_ = FUNCTION_SPACE_HVOL;
206 const ordinal_type card = this->basisCardinality_;
209 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
210 dofCoords(
"HVOL::Line::Cn::dofCoords", card, 1);
214 case POINTTYPE_EQUISPACED:
215 case POINTTYPE_WARPBLEND: {
218 const ordinal_type offset = 1;
220 this->basisCellTopology_,
221 order+1+offset, offset,
240 case POINTTYPE_GAUSS: {
247 INTREPID2_TEST_FOR_EXCEPTION( !isValidPointType(pointType),
248 std::invalid_argument ,
249 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM) invalid pointType." );
253 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoords);
254 Kokkos::deep_copy(this->dofCoords_, dofCoords);
258 const ordinal_type lwork = card*card;
259 Kokkos::DynRankView<typename ScalarViewType::value_type,Kokkos::LayoutLeft,Kokkos::HostSpace>
260 vmat(
"HVOL::Line::Cn::vmat", card, card),
261 work(
"HVOL::Line::Cn::work", lwork),
262 ipiv(
"HVOL::Line::Cn::ipiv", card);
264 const double alpha = 0.0, beta = 0.0;
265 Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
266 getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>
267 (vmat, dofCoords, order, alpha, beta, OPERATOR_VALUE);
269 ordinal_type info = 0;
270 Teuchos::LAPACK<ordinal_type,typename ScalarViewType::value_type> lapack;
272 lapack.GETRF(card, card,
273 vmat.data(), vmat.stride_1(),
274 (ordinal_type*)ipiv.data(),
277 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
279 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM) lapack.GETRF returns nonzero info." );
282 vmat.data(), vmat.stride_1(),
283 (ordinal_type*)ipiv.data(),
287 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
289 ">>> ERROR: (Intrepid2::Basis_HVOL_LINE_Cn_FEM) lapack.GETRI returns nonzero info." );
292 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
293 vinv(
"HVOL::Line::Cn::vinv", card, card);
295 for (ordinal_type i=0;i<card;++i)
296 for (ordinal_type j=0;j<card;++j)
297 vinv(i,j) = vmat(j,i);
299 this->vinv_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), vinv);
300 Kokkos::deep_copy(this->vinv_ , vinv);
305 const ordinal_type tagSize = 4;
306 const ordinal_type posScDim = 0;
307 const ordinal_type posScOrd = 1;
308 const ordinal_type posDfOrd = 2;
313 for (ordinal_type i=0;i<card;++i) {
324 this->setOrdinalTagData(this->tagToOrdinal_,
327 this->basisCardinality_,
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Basis_HVOL_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.