16 #ifndef __INTREPID2_HGRAD_QUAD_C1_FEM_DEF_HPP__
17 #define __INTREPID2_HGRAD_QUAD_C1_FEM_DEF_HPP__
25 template<EOperator opType>
26 template<
typename OutputViewType,
27 typename inputViewType>
28 KOKKOS_INLINE_FUNCTION
30 Basis_HGRAD_QUAD_C1_FEM::Serial<opType>::
31 getValues( OutputViewType output,
32 const inputViewType input ) {
34 case OPERATOR_VALUE : {
35 const auto x = input(0);
36 const auto y = input(1);
39 output.access(0) = (1.0 - x)*(1.0 - y)/4.0;
40 output.access(1) = (1.0 + x)*(1.0 - y)/4.0;
41 output.access(2) = (1.0 + x)*(1.0 + y)/4.0;
42 output.access(3) = (1.0 - x)*(1.0 + y)/4.0;
45 case OPERATOR_GRAD : {
46 const auto x = input(0);
47 const auto y = input(1);
50 output.access(0, 0) = -(1.0 - y)/4.0;
51 output.access(0, 1) = -(1.0 - x)/4.0;
53 output.access(1, 0) = (1.0 - y)/4.0;
54 output.access(1, 1) = -(1.0 + x)/4.0;
56 output.access(2, 0) = (1.0 + y)/4.0;
57 output.access(2, 1) = (1.0 + x)/4.0;
59 output.access(3, 0) = -(1.0 + y)/4.0;
60 output.access(3, 1) = (1.0 - x)/4.0;
63 case OPERATOR_CURL : {
64 const auto x = input(0);
65 const auto y = input(1);
68 output.access(0, 0) = -(1.0 - x)/4.0;
69 output.access(0, 1) = (1.0 - y)/4.0;
71 output.access(1, 0) = -(1.0 + x)/4.0;
72 output.access(1, 1) = -(1.0 - y)/4.0;
74 output.access(2, 0) = (1.0 + x)/4.0;
75 output.access(2, 1) = -(1.0 + y)/4.0;
77 output.access(3, 0) = (1.0 - x)/4.0;
78 output.access(3, 1) = (1.0 + y)/4.0;
83 output.access(0, 0) = 0.0;
84 output.access(0, 1) = 0.25;
85 output.access(0, 2) = 0.0;
87 output.access(1, 0) = 0.0;
88 output.access(1, 1) = -0.25;
89 output.access(1, 2) = 0.0;
91 output.access(2, 0) = 0.0;
92 output.access(2, 1) = 0.25;
93 output.access(2, 2) = 0.0;
95 output.access(3, 0) = 0.0;
96 output.access(3, 1) = -0.25;
97 output.access(3, 2) = 0.0;
100 case OPERATOR_MAX : {
101 const ordinal_type jend = output.extent(1);
102 const ordinal_type iend = output.extent(0);
104 for (ordinal_type j=0;j<jend;++j)
105 for (ordinal_type i=0;i<iend;++i)
106 output.access(i, j) = 0.0;
110 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
111 opType != OPERATOR_GRAD &&
112 opType != OPERATOR_CURL &&
113 opType != OPERATOR_D2 &&
114 opType != OPERATOR_MAX,
115 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::Serial::getValues) operator is not supported");
121 template<
typename DT,
122 typename outputValueValueType,
class ...outputValueProperties,
123 typename inputPointValueType,
class ...inputPointProperties>
125 Basis_HGRAD_QUAD_C1_FEM::
126 getValues(
const typename DT::execution_space& space,
127 Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
128 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
129 const EOperator operatorType ) {
130 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
131 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
132 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
135 const auto loopSize = inputPoints.extent(0);
136 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
138 switch (operatorType) {
140 case OPERATOR_VALUE: {
141 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_VALUE> FunctorType;
142 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
148 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_GRAD> FunctorType;
149 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
152 case OPERATOR_CURL: {
153 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_CURL> FunctorType;
154 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
158 INTREPID2_TEST_FOR_EXCEPTION( operatorType == OPERATOR_DIV, std::invalid_argument,
159 ">>> ERROR (Basis_HGRAD_QUAD_C1_FEM): DIV is invalid operator for rank-0 (scalar) functions in 2D");
163 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_D2> FunctorType;
164 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
175 typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_MAX> FunctorType;
176 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
180 INTREPID2_TEST_FOR_EXCEPTION( !Intrepid2::isValidOperator(operatorType), std::invalid_argument,
181 ">>> ERROR (Basis_HGRAD_QUAD_C1_FEM): Invalid operator type");
188 template<
typename DT,
typename OT,
typename PT>
191 const ordinal_type spaceDim = 2;
192 this->basisCardinality_ = 4;
193 this->basisDegree_ = 1;
194 this->basisCellTopologyKey_ = shards::Quadrilateral<4>::key;
195 this->basisType_ = BASIS_FEM_DEFAULT;
196 this->basisCoordinates_ = COORDINATES_CARTESIAN;
197 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
202 const ordinal_type tagSize = 4;
203 const ordinal_type posScDim = 0;
204 const ordinal_type posScOrd = 1;
205 const ordinal_type posDfOrd = 2;
208 ordinal_type tags[16] = { 0, 0, 0, 1,
219 this->setOrdinalTagData(this->tagToOrdinal_,
222 this->basisCardinality_,
236 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
237 dofCoords(
"dofCoordsHost", this->basisCardinality_,spaceDim);
239 dofCoords(0,0) = -1.0; dofCoords(0,1) = -1.0;
240 dofCoords(1,0) = 1.0; dofCoords(1,1) = -1.0;
241 dofCoords(2,0) = 1.0; dofCoords(2,1) = 1.0;
242 dofCoords(3,0) = -1.0; dofCoords(3,1) = 1.0;
244 this->dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
245 Kokkos::deep_copy(this->dofCoords_, dofCoords);
248 template<
typename DT,
typename OT,
typename PT>
251 ordinal_type& perTeamSpaceSize,
252 ordinal_type& perThreadSpaceSize,
254 const EOperator operatorType)
const {
255 perTeamSpaceSize = 0;
256 perThreadSpaceSize = 0;
259 template<
typename DT,
typename OT,
typename PT>
260 KOKKOS_INLINE_FUNCTION
263 OutputViewType outputValues,
264 const PointViewType inputPoints,
265 const EOperator operatorType,
266 const typename Kokkos::TeamPolicy<typename DT::execution_space>::member_type& team_member,
267 const typename DT::execution_space::scratch_memory_space & scratchStorage,
268 const ordinal_type subcellDim,
269 const ordinal_type subcellOrdinal)
const {
271 INTREPID2_TEST_FOR_ABORT( !((subcellDim <= 0) && (subcellOrdinal == -1)),
272 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getValues), The capability of selecting subsets of basis functions has not been implemented yet.");
274 (void) scratchStorage;
276 const int numPoints = inputPoints.extent(0);
278 switch(operatorType) {
280 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
281 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
282 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
287 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
288 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
289 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
290 Impl::Basis_HGRAD_QUAD_C1_FEM::Serial<OPERATOR_GRAD>::getValues( output, input);
294 Kokkos::parallel_for (Kokkos::TeamThreadRange (team_member, numPoints), [=] (ordinal_type& pt) {
295 auto output = Kokkos::subview( outputValues, Kokkos::ALL(), pt, Kokkos::ALL() );
296 const auto input = Kokkos::subview( inputPoints, pt, Kokkos::ALL() );
297 Impl::Basis_HGRAD_QUAD_C1_FEM::Serial<OPERATOR_CURL>::getValues( output, input);
301 INTREPID2_TEST_FOR_ABORT(
true,
">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::getValues), Operator Type not supported.");
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...
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.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Kokkos::DynRankView< PointValueType, Kokkos::LayoutStride, DeviceType > PointViewType
View type for input points.
Basis_HGRAD_QUAD_C1_FEM()
Constructor.
See Intrepid2::Basis_HGRAD_QUAD_C1_FEM.