16 #ifndef __INTREPID2_HGRAD_TRI_C1_FEM_DEF_HPP__
17 #define __INTREPID2_HGRAD_TRI_C1_FEM_DEF_HPP__
26 template<EOperator opType>
27 template<
typename OutputViewType,
28 typename inputViewType>
29 KOKKOS_INLINE_FUNCTION
31 Basis_HGRAD_TRI_C1_FEM::Serial<opType>::
32 getValues( OutputViewType output,
33 const inputViewType input ) {
35 case OPERATOR_VALUE: {
36 const auto x = input(0);
37 const auto y = input(1);
40 output.access(0) = 1.0 - x - y;
47 output.access(0, 0) = -1.0;
48 output.access(0, 1) = -1.0;
50 output.access(1, 0) = 1.0;
51 output.access(1, 1) = 0.0;
53 output.access(2, 0) = 0.0;
54 output.access(2, 1) = 1.0;
58 output.access(0, 0) = -1.0;
59 output.access(0, 1) = 1.0;
61 output.access(1, 0) = 0.0;
62 output.access(1, 1) = -1.0;
64 output.access(2, 0) = 1.0;
65 output.access(2, 1) = 0.0;
69 const ordinal_type jend = output.extent(1);
70 const ordinal_type iend = output.extent(0);
72 for (ordinal_type j=0;j<jend;++j)
73 for (ordinal_type i=0;i<iend;++i)
74 output.access(i, j) = 0.0;
78 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
79 opType != OPERATOR_GRAD &&
80 opType != OPERATOR_CURL &&
81 opType != OPERATOR_MAX,
82 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::Serial::getValues) operator is not supported");
89 typename outputValueValueType,
class ...outputValueProperties,
90 typename inputPointValueType,
class ...inputPointProperties>
92 Basis_HGRAD_TRI_C1_FEM::
93 getValues(
const typename DT::execution_space& space,
94 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
95 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
96 const EOperator operatorType ) {
97 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
98 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
99 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
102 const auto loopSize = inputPoints.extent(0);
103 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
105 switch (operatorType) {
106 case OPERATOR_VALUE: {
107 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_VALUE> FunctorType;
108 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
113 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_GRAD> FunctorType;
114 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
117 case OPERATOR_CURL: {
118 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_CURL> FunctorType;
119 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
123 INTREPID2_TEST_FOR_EXCEPTION( operatorType == OPERATOR_DIV, std::invalid_argument,
124 ">>> ERROR (Basis_HGRAD_TRI_C1_FEM): DIV is invalid operator for rank-0 (scalar) functions in 2D");
136 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_MAX> FunctorType;
137 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
141 INTREPID2_TEST_FOR_EXCEPTION( !Intrepid2::isValidOperator(operatorType), std::invalid_argument,
142 ">>> ERROR (Basis_HGRAD_TRI_C1_FEM): Invalid operator type");
149 template<
typename DT,
typename OT,
typename PT>
152 const ordinal_type spaceDim = 2;
153 this->basisCardinality_ = 3;
154 this->basisDegree_ = 1;
155 this->basisCellTopologyKey_ = shards::Triangle<3>::key;
156 this->basisType_ = BASIS_FEM_DEFAULT;
157 this->basisCoordinates_ = COORDINATES_CARTESIAN;
158 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
163 const ordinal_type tagSize = 4;
164 const ordinal_type posScDim = 0;
165 const ordinal_type posScOrd = 1;
166 const ordinal_type posDfOrd = 2;
169 ordinal_type tags[12] = { 0, 0, 0, 1,
179 this->setOrdinalTagData(this->tagToOrdinal_,
182 this->basisCardinality_,
196 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
197 dofCoords(
"dofCoordsHost", this->basisCardinality_,spaceDim);
199 dofCoords(0,0) = 0.0; dofCoords(0,1) = 0.0;
200 dofCoords(1,0) = 1.0; dofCoords(1,1) = 0.0;
201 dofCoords(2,0) = 0.0; dofCoords(2,1) = 1.0;
203 this->dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
204 Kokkos::deep_copy(this->dofCoords_, dofCoords);
207 template<
typename DT,
typename OT,
typename PT>
210 ordinal_type& perTeamSpaceSize,
211 ordinal_type& perThreadSpaceSize,
213 const EOperator operatorType)
const {
214 perTeamSpaceSize = 0;
215 perThreadSpaceSize = 0;
218 template<
typename DT,
typename OT,
typename PT>
219 KOKKOS_INLINE_FUNCTION
222 OutputViewType outputValues,
223 const PointViewType inputPoints,
224 const EOperator operatorType,
225 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
226 const typename DT::execution_space::scratch_memory_space & scratchStorage,
227 const ordinal_type subcellDim,
228 const ordinal_type subcellOrdinal)
const {
230 INTREPID2_TEST_FOR_ABORT( !((subcellDim <= 0) && (subcellOrdinal == -1)),
231 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
233 (void) scratchStorage;
235 const int numPoints = inputPoints.extent(0);
237 switch(operatorType) {
239 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
240 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
241 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
246 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
247 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
248 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
249 Impl::Basis_HGRAD_TRI_C1_FEM::Serial<OPERATOR_GRAD>::getValues( output, input);
253 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
254 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
255 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
256 Impl::Basis_HGRAD_TRI_C1_FEM::Serial<OPERATOR_CURL>::getValues( output, input);
260 INTREPID2_TEST_FOR_ABORT(
true,
">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_C1_FEM::getValues), Operator Type not supported.");
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
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.
See Intrepid2::Basis_HGRAD_TRI_C1_FEM.
Basis_HGRAD_TRI_C1_FEM()
Constructor.
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...