Intrepid2
Intrepid2_HGRAD_LINE_C2_FEMDef.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Intrepid2 Package
4//
5// Copyright 2007 NTESS and the Intrepid2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14#ifndef __INTREPID2_HGRAD_LINE_C2_FEM_DEF_HPP__
15#define __INTREPID2_HGRAD_LINE_C2_FEM_DEF_HPP__
16
17namespace Intrepid2 {
18
19 // -------------------------------------------------------------------------------------
20
21 namespace Impl {
22
23 template<EOperator opType>
24 template<typename OutputViewType,
25 typename inputViewType>
26 KOKKOS_INLINE_FUNCTION
27 void
28 Basis_HGRAD_LINE_C2_FEM::Serial<opType>::
29 getValues( OutputViewType output,
30 const inputViewType input ) {
31 switch (opType) {
32 case OPERATOR_VALUE : {
33 const auto x = input(0);
34
35 output.access(0) = (x - 1.0)*x/2.0;
36 output.access(1) = (1.0 + x)*x/2.0;
37 output.access(2) = (1.0 + x)*(1.0 - x);
38 break;
39 }
40 case OPERATOR_GRAD : {
41 const auto x = input(0);
42
43 output.access(0, 0) = x-0.5;
44 output.access(1, 0) = x+0.5;
45 output.access(2, 0) = -2.0*x;
46 break;
47 }
48 case OPERATOR_MAX : {
49 const ordinal_type jend = output.extent(1);
50 const ordinal_type iend = output.extent(0);
51
52 for (ordinal_type j=0;j<jend;++j)
53 for (ordinal_type i=0;i<iend;++i)
54 output.access(i, j) = 0.0;
55 break;
56 }
57 default: {
58 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
59 opType != OPERATOR_GRAD &&
60 opType != OPERATOR_MAX,
61 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::Serial::getValues) operator is not supported");
62
63 }
64 }
65 }
66
67 template<typename DT,
68 typename outputValueValueType, class ...outputValueProperties,
69 typename inputPointValueType, class ...inputPointProperties>
70 void
71 Basis_HGRAD_LINE_C2_FEM::
72 getValues( const typename DT::execution_space& space,
73 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
74 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
75 const EOperator operatorType ) {
76 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
77 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
78 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
79
80 // Number of evaluation points = dim 0 of inputPoints
81 const auto loopSize = inputPoints.extent(0);
82 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
83
84 switch (operatorType) {
85
86 case OPERATOR_VALUE: {
87 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_VALUE> FunctorType;
88 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
89 break;
90 }
91 case OPERATOR_GRAD:
92 case OPERATOR_DIV:
93 case OPERATOR_CURL:
94 case OPERATOR_D1: {
95 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_GRAD> FunctorType;
96 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
97 break;
98 }
99 case OPERATOR_D2:
100 case OPERATOR_D3:
101 case OPERATOR_D4:
102 case OPERATOR_D5:
103 case OPERATOR_D6:
104 case OPERATOR_D7:
105 case OPERATOR_D8:
106 case OPERATOR_D9:
107 case OPERATOR_D10: {
108 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_MAX> FunctorType;
109 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
110 break;
111 }
112 default: {
113 INTREPID2_TEST_FOR_EXCEPTION( !Intrepid2::isValidOperator(operatorType), std::invalid_argument,
114 ">>> ERROR (Basis_HGRAD_LINE_C2_FEM): Invalid operator type");
115 }
116 }
117 }
118
119
120
121 }
122
123 // -------------------------------------------------------------------------------------
124
125 template<typename DT, typename OT, typename PT>
128 const ordinal_type spaceDim = 1;
129 this->basisCardinality_ = 3;
130 this->basisDegree_ = 2;
131 this->basisCellTopologyKey_ = shards::Line<2>::key;
132 this->basisType_ = BASIS_FEM_DEFAULT;
133 this->basisCoordinates_ = COORDINATES_CARTESIAN;
134 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
135
136 // initialize tags
137 {
138 // Basis-dependent intializations
139 const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
140 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
141 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
142 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
143
144 // An array with local DoF tags assigned to basis functions, in the order of their local enumeration
145 ordinal_type tags[12] = { 0, 0, 0, 1,
146 0, 1, 0, 1,
147 1, 0, 0, 1 };
148
149 // host tags
150 OrdinalTypeArray1DHost tagView(&tags[0], 12);
151
152 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
153 // tags are constructed on host and sent to devices
154 this->setOrdinalTagData(this->tagToOrdinal_,
155 this->ordinalToTag_,
156 tagView,
157 this->basisCardinality_,
158 tagSize,
159 posScDim,
160 posScOrd,
161 posDfOrd);
162 }
163
164 // dofCoords on host and create its mirror view to device
165 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
166 dofCoords("dofCoordsHost", this->basisCardinality_,spaceDim);
167
168 dofCoords(0,0) = -1.0;
169 dofCoords(1,0) = 1.0;
170 dofCoords(2,0) = 0.0;
171
172 this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
173 Kokkos::deep_copy(this->dofCoords_, dofCoords);
174 }
175
176 template<typename DT, typename OT, typename PT>
177 void
179 ordinal_type& perTeamSpaceSize,
180 ordinal_type& perThreadSpaceSize,
181 const PointViewType inputPoints,
182 const EOperator operatorType) const {
183 perTeamSpaceSize = 0;
184 perThreadSpaceSize = 0;
185 }
186
187 template<typename DT, typename OT, typename PT>
188 KOKKOS_INLINE_FUNCTION
189 void
191 OutputViewType outputValues,
192 const PointViewType inputPoints,
193 const EOperator operatorType,
194 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
195 const typename DT::execution_space::scratch_memory_space & scratchStorage,
196 const ordinal_type subcellDim,
197 const ordinal_type subcellOrdinal) const {
198
199 INTREPID2_TEST_FOR_ABORT( !((subcellDim == -1) && (subcellOrdinal == -1)),
200 ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_C2_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
201
202 (void) scratchStorage; //avoid unused variable warning
203
204 const int numPoints = inputPoints.extent(0);
205
206 switch(operatorType) {
207 case OPERATOR_VALUE:
208 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
209 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
210 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
212 });
213 break;
214 case OPERATOR_GRAD:
215 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
216 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
217 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
218 Impl::Basis_HGRAD_LINE_C2_FEM::Serial<OPERATOR_GRAD>::getValues( output, input);
219 });
220 break;
221 default: {}
222 }
223 }
224
225}// namespace Intrepid2
226
227#endif
virtual void getValues(const ExecutionSpace &space, OutputViewType outputValues, const PointViewType inputPoints, const EOperator operatorType=OPERATOR_VALUE) const override
Evaluation of a FEM basis on a reference cell.
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...
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.