16 #ifndef __INTREPID2_HGRAD_TET_CN_FEM_ORTH_DEF_HPP__
17 #define __INTREPID2_HGRAD_TET_CN_FEM_ORTH_DEF_HPP__
24 template<
typename OutputViewType,
25 typename inputViewType,
26 typename workViewType,
28 KOKKOS_INLINE_FUNCTION
29 void OrthPolynomialTet<OutputViewType,inputViewType,workViewType,hasDeriv,0>::generate(
30 OutputViewType output,
31 const inputViewType input,
33 const ordinal_type order ) {
35 constexpr ordinal_type spaceDim = 3;
38 typedef typename OutputViewType::value_type value_type;
40 auto output0 = (hasDeriv) ? Kokkos::subview(output, Kokkos::ALL(), Kokkos::ALL(),0) : Kokkos::subview(output, Kokkos::ALL(), Kokkos::ALL());
43 npts = input.extent(0);
53 const ordinal_type loc = Intrepid2::getPnEnumeration<spaceDim>(0,0,0);
54 for (ordinal_type i=0;i<npts;++i) {
55 output0(loc, i) = 1.0;
57 output.access(loc,i,1) = 0;
58 output.access(loc,i,2) = 0;
59 output.access(loc,i,3) = 0;
65 value_type f1[maxNumPts]={},f2[maxNumPts]={}, f3[maxNumPts]={}, f4[maxNumPts]={},f5[maxNumPts]={},
66 df2_1[maxNumPts]={}, df2_2[maxNumPts]={}, df5_2[maxNumPts]={};
67 value_type df1_0, df1_1, df1_2, df3_1, df3_2, df4_2;
69 for (
int i=0;i<npts;i++) {
70 f1[i] = 2.0*z(i,0) + z(i,1) + z(i,2)- 1.0;
71 f2[i] = pow(z(i,1) + z(i,2) -1.0, 2);
72 f3[i] = 2.0*z(i,1) + z(i,2) -1.0;
74 f5[i] = f4[i] * f4[i];
76 df2_1[i] = 2*(z(i,1) + z(i,2) -1.0);
91 const ordinal_type loc = Intrepid2::getPnEnumeration<spaceDim>(1,0,0);
92 for (ordinal_type i=0;i<npts;++i) {
93 output0(loc, i) = f1[i];
95 output.access(loc,i,1) = df1_0;
96 output.access(loc,i,2) = df1_1;
97 output.access(loc,i,3) = df1_2;
104 for (ordinal_type p=1;p<order;p++) {
106 loc = Intrepid2::getPnEnumeration<spaceDim>(p,0,0),
107 loc_p1 = Intrepid2::getPnEnumeration<spaceDim>(p+1,0,0),
108 loc_m1 = Intrepid2::getPnEnumeration<spaceDim>(p-1,0,0);
111 a = (2.0*p+1.0)/(1.0+p),
114 for (ordinal_type i=0;i<npts;++i) {
115 output0(loc_p1,i) = ( a * f1[i] * output0(loc,i) -
116 b * f2[i] * output0(loc_m1,i) );
118 output.access(loc_p1,i,1) = a * (f1[i] * output.access(loc,i,1) + df1_0 * output0(loc,i)) -
119 b * f2[i] * output.access(loc_m1,i,1) ;
120 output.access(loc_p1,i,2) = a * (f1[i] * output.access(loc,i,2) + df1_1 * output0(loc,i)) -
121 b * (df2_1[i] * output0(loc_m1,i) + f2[i] * output.access(loc_m1,i,2)) ;
122 output.access(loc_p1,i,3) = a * (f1[i] * output.access(loc,i,3) + df1_2 * output0(loc,i)) -
123 b * (df2_2[i] * output0(loc_m1,i) + f2[i] * output.access(loc_m1,i,3)) ;
129 for (ordinal_type p=0;p<order;++p) {
131 loc_p_0 = Intrepid2::getPnEnumeration<spaceDim>(p,0,0),
132 loc_p_1 = Intrepid2::getPnEnumeration<spaceDim>(p,1,0);
134 for (ordinal_type i=0;i<npts;++i) {
135 const value_type coeff = (2.0 * p + 3.0) * z(i,1) + z(i,2)-1.0;
136 output0(loc_p_1,i) = output0(loc_p_0,i)*coeff;
138 output.access(loc_p_1,i,1) = output.access(loc_p_0,i,1)*coeff;
139 output.access(loc_p_1,i,2) = output.access(loc_p_0,i,2)*coeff + output0(loc_p_0,i)*(2.0 * p + 3.0);
140 output.access(loc_p_1,i,3) = output.access(loc_p_0,i,3)*coeff + output0(loc_p_0,i);
148 for (ordinal_type p=0;p<order-1;++p)
149 for (ordinal_type q=1;q<order-p;++q) {
151 loc_p_qp1 = Intrepid2::getPnEnumeration<spaceDim>(p,q+1,0),
152 loc_p_q = Intrepid2::getPnEnumeration<spaceDim>(p,q,0),
153 loc_p_qm1 = Intrepid2::getPnEnumeration<spaceDim>(p,q-1,0);
155 value_type a,b,c, coeff, dcoeff_1, dcoeff_2;
156 Intrepid2::getJacobyRecurrenceCoeffs(a,b,c, 2*p+1,0,q);
157 for (ordinal_type i=0;i<npts;++i) {
158 coeff = a * f3[i] + b * f4[i];
159 dcoeff_1 = a * df3_1;
160 dcoeff_2 = a * df3_2 + b * df4_2;
161 output0(loc_p_qp1,i) = coeff * output0(loc_p_q,i)
162 - c* f5[i] * output0(loc_p_qm1,i) ;
164 output.access(loc_p_qp1,i,1) = coeff * output.access(loc_p_q,i,1) +
165 - c * f5[i] * output.access(loc_p_qm1,i,1) ;
166 output.access(loc_p_qp1,i,2) = coeff * output.access(loc_p_q,i,2) + dcoeff_1 * output0(loc_p_q,i)
167 - c * f5[i] * output.access(loc_p_qm1,i,2) ;
168 output.access(loc_p_qp1,i,3) = coeff * output.access(loc_p_q,i,3) + dcoeff_2 * output0(loc_p_q,i)
169 - c * f5[i] * output.access(loc_p_qm1,i,3) - c * df5_2[i] * output0(loc_p_qm1,i);
176 for (ordinal_type p=0;p<order;++p)
177 for (ordinal_type q=0;q<order-p;++q) {
180 loc_p_q_0 = Intrepid2::getPnEnumeration<spaceDim>(p,q,0),
181 loc_p_q_1 = Intrepid2::getPnEnumeration<spaceDim>(p,q,1);
183 for (ordinal_type i=0;i<npts;++i) {
184 const value_type coeff = 2.0 * ( 2.0 + q + p ) * z(i,2) - 1.0;
185 output0(loc_p_q_1,i) = output0(loc_p_q_0,i) * coeff;
187 output.access(loc_p_q_1,i,1) = output.access(loc_p_q_0,i,1) * coeff;
188 output.access(loc_p_q_1,i,2) = output.access(loc_p_q_0,i,2) * coeff;
189 output.access(loc_p_q_1,i,3) = output.access(loc_p_q_0,i,3) * coeff + 2.0 * ( 2.0 + q + p ) * output0(loc_p_q_0,i);
196 for (ordinal_type p=0;p<order-1;++p)
197 for (ordinal_type q=0;q<order-p-1;++q)
198 for (ordinal_type r=1;r<order-p-q;++r) {
200 loc_p_q_rp1 = Intrepid2::getPnEnumeration<spaceDim>(p,q,r+1),
201 loc_p_q_r = Intrepid2::getPnEnumeration<spaceDim>(p,q,r),
202 loc_p_q_rm1 = Intrepid2::getPnEnumeration<spaceDim>(p,q,r-1);
204 value_type a,b,c, coeff;
205 Intrepid2::getJacobyRecurrenceCoeffs(a,b,c, 2*p+2*q+2,0,r);
206 for (ordinal_type i=0;i<npts;++i) {
207 coeff = 2.0 * a * z(i,2) - a + b;
208 output0(loc_p_q_rp1,i) = coeff * output0(loc_p_q_r,i) - c * output0(loc_p_q_rm1,i) ;
210 output.access(loc_p_q_rp1,i,1) = coeff * output.access(loc_p_q_r,i,1) - c * output.access(loc_p_q_rm1,i,1);
211 output.access(loc_p_q_rp1,i,2) = coeff * output.access(loc_p_q_r,i,2) - c * output.access(loc_p_q_rm1,i,2);
212 output.access(loc_p_q_rp1,i,3) = coeff * output.access(loc_p_q_r,i,3) + 2 * a * output0(loc_p_q_r,i) - c * output.access(loc_p_q_rm1,i,3);
220 for (ordinal_type p=0;p<=order;++p)
221 for (ordinal_type q=0;q<=order-p;++q)
222 for (ordinal_type r=0;r<=order-p-q;++r) {
223 ordinal_type loc_p_q_r = Intrepid2::getPnEnumeration<spaceDim>(p,q,r);
224 value_type scal = std::sqrt( (p+0.5)*(p+q+1.0)*(p+q+r+1.5) );
225 for (ordinal_type i=0;i<npts;++i) {
226 output0(loc_p_q_r,i) *= scal;
228 output.access(loc_p_q_r,i,1) *= scal;
229 output.access(loc_p_q_r,i,2) *= scal;
230 output.access(loc_p_q_r,i,3) *= scal;
236 template<
typename OutputViewType,
237 typename inputViewType,
238 typename workViewType,
240 KOKKOS_INLINE_FUNCTION
241 void OrthPolynomialTet<OutputViewType,inputViewType,workViewType,hasDeriv,1>::generate(
242 OutputViewType output,
243 const inputViewType input,
245 const ordinal_type order ) {
246 constexpr ordinal_type spaceDim = 3;
248 npts = input.extent(0),
249 card = output.extent(0);
251 workViewType dummyView;
252 OrthPolynomialTet<workViewType,inputViewType,workViewType,hasDeriv,0>::generate(work, input, dummyView, order);
253 for (ordinal_type i=0;i<card;++i)
254 for (ordinal_type j=0;j<npts;++j)
255 for (ordinal_type k=0;k<spaceDim;++k)
256 output.access(i,j,k) = work(i,j,k+1);
261 template<
typename OutputViewType,
262 typename inputViewType,
263 typename workViewType,
266 KOKKOS_INLINE_FUNCTION
267 void OrthPolynomialTet<OutputViewType,inputViewType,workViewType,hasDeriv,n>::generate(
269 const inputViewType ,
271 const ordinal_type ) {
272 #if 0 //#ifdef HAVE_INTREPID2_SACADO
274 constexpr ordinal_type spaceDim = 3;
275 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
277 typedef typename OutputViewType::value_type value_type;
278 typedef Sacado::Fad::SFad<value_type,spaceDim> fad_type;
281 npts = input.extent(0),
282 card = output.extent(0);
288 typedef typename inputViewType::memory_space memory_space;
289 typedef typename inputViewType::memory_space memory_space;
290 typedef typename Kokkos::View<fad_type***, memory_space> outViewType;
291 typedef typename Kokkos::View<fad_type**, memory_space> inViewType;
292 auto vcprop = Kokkos::common_view_alloc_prop(input);
294 inViewType in(Kokkos::view_wrap((value_type*)&inBuf[0][0], vcprop), npts, spaceDim);
295 outViewType out(Kokkos::view_wrap((value_type*)&outBuf[0][0][0], vcprop), card, npts, n*(n+1)/2);
297 for (ordinal_type i=0;i<npts;++i)
298 for (ordinal_type j=0;j<spaceDim;++j) {
299 in(i,j) = input(i,j);
300 in(i,j).diff(j,spaceDim);
303 typedef typename Kokkos::DynRankView<fad_type, memory_space> outViewType_;
304 outViewType_ workView;
308 auto vcprop = Kokkos::common_view_alloc_prop(in);
309 workView = outViewType_( Kokkos::view_wrap((value_type*)&outBuf[0][0][0], vcprop), card, npts, spaceDim+1);
311 OrthPolynomialTet<outViewType,inViewType,outViewType_,hasDeriv,n-1>::generate(out, in, workView, order);
313 for (ordinal_type i=0;i<card;++i)
314 for (ordinal_type j=0;j<npts;++j) {
315 for (ordinal_type i_dx = 0; i_dx <= n; ++i_dx)
316 for (ordinal_type i_dy = 0; i_dy <= n-i_dx; ++i_dy) {
317 ordinal_type i_dz = n - i_dx - i_dy;
318 ordinal_type i_Dn = Intrepid2::getDkEnumeration<spaceDim>(i_dx,i_dy,i_dz);
322 ordinal_type i_Dnm1 = Intrepid2::getDkEnumeration<spaceDim>(i_dx-1, i_dy, i_dz);
323 output.access(i,j,i_Dn) = out(i,j,i_Dnm1).dx(0);
328 ordinal_type i_Dnm1 = Intrepid2::getDkEnumeration<spaceDim>(i_dx, i_dy-1, i_dz);
329 output.access(i,j,i_Dn) = out(i,j,i_Dnm1).dx(1);
334 ordinal_type i_Dnm1 = Intrepid2::getDkEnumeration<spaceDim>(i_dx, i_dy, i_dz-1);
335 output.access(i,j,i_Dn) = out(i,j,i_Dnm1).dx(2);
340 INTREPID2_TEST_FOR_ABORT(
true,
341 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH::OrthPolynomialTri) Computing of second and higher-order derivatives is not currently supported");
346 template<EOperator opType>
347 template<
typename OutputViewType,
348 typename inputViewType,
349 typename workViewType>
350 KOKKOS_INLINE_FUNCTION
352 Basis_HGRAD_TET_Cn_FEM_ORTH::Serial<opType>::
353 getValues( OutputViewType output,
354 const inputViewType input,
356 const ordinal_type order) {
358 case OPERATOR_VALUE: {
359 OrthPolynomialTet<OutputViewType,inputViewType,workViewType,false,0>::generate( output, input, work, order );
364 OrthPolynomialTet<OutputViewType,inputViewType,workViewType,true,1>::generate( output, input, work, order );
368 OrthPolynomialTet<OutputViewType,inputViewType,workViewType,true,2>::generate( output, input, work, order );
372 INTREPID2_TEST_FOR_ABORT(
true,
373 ">>> ERROR: (Intrepid2::Basis_HGRAD_TET_Cn_FEM_ORTH::Serial::getValues) operator is not supported");
380 template<
typename DT, ordinal_type numPtsPerEval,
381 typename outputValueValueType,
class ...outputValueProperties,
382 typename inputPointValueType,
class ...inputPointProperties>
384 Basis_HGRAD_TET_Cn_FEM_ORTH::
386 const typename DT::execution_space& space,
387 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
388 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
389 const ordinal_type order,
390 const EOperator operatorType ) {
391 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
392 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
393 typedef typename DT::execution_space ExecSpaceType;
396 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
397 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
398 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
399 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
401 typedef typename inputPointViewType::value_type inputPointType;
402 const ordinal_type cardinality = outputValues.extent(0);
403 const ordinal_type spaceDim = 3;
405 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
406 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
408 switch (operatorType) {
409 case OPERATOR_VALUE: {
410 workViewType dummyWorkView;
411 typedef Functor<outputValueViewType,inputPointViewType,workViewType,OPERATOR_VALUE,numPtsPerEval> FunctorType;
412 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, dummyWorkView, order) );
417 workViewType work(Kokkos::view_alloc(space,
"Basis_HGRAD_TET_In_FEM_ORTH::getValues::work", vcprop), cardinality, inputPoints.extent(0), spaceDim+1);
418 typedef Functor<outputValueViewType,inputPointViewType,workViewType,OPERATOR_D1,numPtsPerEval> FunctorType;
419 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, work, order) );
423 workViewType dummyWorkView;
424 typedef Functor<outputValueViewType,inputPointViewType,workViewType,OPERATOR_D2,numPtsPerEval> FunctorType;
425 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints ,dummyWorkView, order) );
429 INTREPID2_TEST_FOR_EXCEPTION( !Intrepid2::isValidOperator(operatorType), std::invalid_argument,
430 ">>> ERROR (Basis_HGRAD_TET_Cn_FEM_ORTH): invalid operator type");
438 template<
typename DT,
typename OT,
typename PT>
442 constexpr ordinal_type spaceDim = 3;
443 this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order);
444 this->basisDegree_ = order;
445 this->basisCellTopologyKey_ = shards::Tetrahedron<4>::key;
446 this->basisType_ = BASIS_FEM_HIERARCHICAL;
447 this->basisCoordinates_ = COORDINATES_CARTESIAN;
448 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
453 constexpr ordinal_type tagSize = 4;
454 const ordinal_type posScDim = 0;
455 const ordinal_type posScOrd = 1;
456 const ordinal_type posDfOrd = 2;
458 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
459 ordinal_type tags[maxCard][tagSize];
460 const ordinal_type card = this->basisCardinality_;
461 for (ordinal_type i=0;i<card;++i) {
472 this->setOrdinalTagData(this->tagToOrdinal_,
475 this->basisCardinality_,
static constexpr ordinal_type MaxNumPtsPerBasisEval
The maximum number of points to eval in serial mode.
Basis_HGRAD_TET_Cn_FEM_ORTH(const ordinal_type order)
Constructor.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.