49 #ifndef __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
50 #define __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
59 template<EOperator opType>
60 template<
typename OutputViewType,
61 typename inputViewType,
62 typename workViewType,
63 typename vinvViewType>
64 KOKKOS_INLINE_FUNCTION
66 Basis_HGRAD_TRI_Cn_FEM::Serial<opType>::
67 getValues( OutputViewType output,
68 const inputViewType input,
70 const vinvViewType vinv ) {
72 constexpr ordinal_type spaceDim = 2;
74 card = vinv.extent(0),
75 npts = input.extent(0);
78 ordinal_type order = 0;
80 if (card == Intrepid2::getPnCardinality<spaceDim>(p) ) {
86 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
87 auto vcprop = Kokkos::common_view_alloc_prop(work);
88 auto ptr = work.data();
91 case OPERATOR_VALUE: {
92 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
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 output.access(i,j) = 0.0;
101 for (ordinal_type k=0;k<card;++k)
102 output.access(i,j) += vinv(k,i)*phis.access(k,j);
108 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
109 ptr += card*npts*spaceDim*get_dimension_scalar(work);
110 const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
111 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
112 Serial<opType>::getValues(phis, input, workView, order);
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<spaceDim;++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);
132 const ordinal_type dkcard = getDkCardinality<opType,spaceDim>();
133 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
136 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
137 Serial<opType>::getValues(phis, input, dummyView, order);
139 for (ordinal_type i=0;i<card;++i)
140 for (ordinal_type j=0;j<npts;++j)
141 for (ordinal_type k=0;k<dkcard;++k) {
142 output.access(i,j,k) = 0.0;
143 for (ordinal_type l=0;l<card;++l)
144 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
148 case OPERATOR_CURL: {
149 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
150 ptr += card*npts*spaceDim*get_dimension_scalar(work);
151 const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
154 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
155 Serial<OPERATOR_D1>::getValues(phis, input, workView, order);
157 for (ordinal_type i=0;i<card;++i)
158 for (ordinal_type j=0;j<npts;++j) {
159 output.access(i,j,0) = 0.0;
160 for (ordinal_type l=0;l<card;++l)
161 output.access(i,j,0) += vinv(l,i)*phis.access(l,j,1);
162 output.access(i,j,1) = 0.0;
163 for (ordinal_type l=0;l<card;++l)
164 output.access(i,j,1) -= vinv(l,i)*phis.access(l,j,0);
169 INTREPID2_TEST_FOR_ABORT(
true,
170 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented");
175 template<
typename SpT, ordinal_type numPtsPerEval,
176 typename outputValueValueType,
class ...outputValueProperties,
177 typename inputPointValueType,
class ...inputPointProperties,
178 typename vinvValueType,
class ...vinvProperties>
180 Basis_HGRAD_TRI_Cn_FEM::
181 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
182 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
183 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
184 const EOperator operatorType) {
185 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
186 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
187 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
188 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
191 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
192 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
193 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
194 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
196 typedef typename inputPointViewType::value_type inputPointType;
198 const ordinal_type cardinality = outputValues.extent(0);
199 const ordinal_type spaceDim = 2;
201 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
202 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
204 switch (operatorType) {
205 case OPERATOR_VALUE: {
206 workViewType work(Kokkos::view_alloc(
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
207 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
208 OPERATOR_VALUE,numPtsPerEval> FunctorType;
209 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
214 workViewType work(Kokkos::view_alloc(
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
215 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
216 OPERATOR_D1,numPtsPerEval> FunctorType;
217 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
220 case OPERATOR_CURL: {
221 workViewType work(Kokkos::view_alloc(
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
222 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
223 OPERATOR_CURL,numPtsPerEval> FunctorType;
224 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
228 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
229 OPERATOR_D2,numPtsPerEval> FunctorType;
230 workViewType work(Kokkos::view_alloc(
"Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
231 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
242 INTREPID2_TEST_FOR_EXCEPTION(
true , std::invalid_argument,
243 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented" );
250 template<
typename SpT,
typename OT,
typename PT>
253 const EPointType pointType ) {
254 constexpr ordinal_type spaceDim = 2;
256 this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order);
257 this->basisDegree_ = order;
258 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
259 this->basisType_ = BASIS_FEM_FIAT;
260 this->basisCoordinates_ = COORDINATES_CARTESIAN;
261 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
263 const ordinal_type card = this->basisCardinality_;
266 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
267 dofCoords(
"Hgrad::Tri::Cn::dofCoords", card, spaceDim);
270 const ordinal_type offset = 0;
272 this->basisCellTopology_,
276 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoords);
277 Kokkos::deep_copy(this->dofCoords_, dofCoords);
281 const ordinal_type lwork = card*card;
282 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
283 vmat(
"Hgrad::Tri::Cn::vmat", card, card),
284 work(
"Hgrad::Tri::Cn::work", lwork),
285 ipiv(
"Hgrad::Tri::Cn::ipiv", card);
287 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(vmat, dofCoords, order, OPERATOR_VALUE);
289 ordinal_type info = 0;
290 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
292 lapack.GETRF(card, card,
293 vmat.data(), vmat.stride_1(),
294 (ordinal_type*)ipiv.data(),
297 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
299 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
302 vmat.data(), vmat.stride_1(),
303 (ordinal_type*)ipiv.data(),
307 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
309 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
312 Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
313 vinv(
"Hgrad::Line::Cn::vinv", card, card);
315 for (ordinal_type i=0;i<card;++i)
316 for (ordinal_type j=0;j<card;++j)
317 vinv(i,j) = vmat(j,i);
319 this->vinv_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), vinv);
320 Kokkos::deep_copy(this->vinv_ , vinv);
325 constexpr ordinal_type tagSize = 4;
326 const ordinal_type posScDim = 0;
327 const ordinal_type posScOrd = 1;
328 const ordinal_type posDfOrd = 2;
330 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
331 ordinal_type tags[maxCard][tagSize];
334 numEdgeDof = Intrepid2::getPnCardinality<1>(order-2),
335 numElemDof = (order > 2 ? Intrepid2::getPnCardinality<2>(order-3) : 0);
337 scalarType xi0, xi1, xi2;
338 const scalarType eps = threshold();
340 ordinal_type edgeId[3] = {}, elemId = 0;
341 for (ordinal_type i=0;i<card;++i) {
344 const auto x = dofCoords(i,0);
345 const auto y = dofCoords(i,1);
351 if ((1.0 - xi0) < eps) {
357 else if ((1.0 - xi1) < eps) {
363 else if ((1.0 - xi2) < eps) {
369 else if (xi2 < eps) {
372 tags[i][2] = edgeId[0]++;
373 tags[i][3] = numEdgeDof;
375 else if (xi0 < eps) {
378 tags[i][2] = edgeId[1]++;
379 tags[i][3] = numEdgeDof;
381 else if (xi1 < eps) {
384 tags[i][2] = edgeId[2]++;
385 tags[i][3] = numEdgeDof;
390 tags[i][2] = elemId++;
391 tags[i][3] = numElemDof;
399 this->setOrdinalTagData(this->tagToOrdinal_,
402 this->basisCardinality_,
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Basis_HGRAD_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.