49 #ifndef __INTREPID2_HCURL_QUAD_I1_FEM_DEF_HPP__
50 #define __INTREPID2_HCURL_QUAD_I1_FEM_DEF_HPP__
58 template<EOperator opType>
59 template<
typename OutputViewType,
60 typename inputViewType>
61 KOKKOS_INLINE_FUNCTION
63 Basis_HCURL_QUAD_I1_FEM::Serial<opType>::
64 getValues( OutputViewType output,
65 const inputViewType input ) {
67 case OPERATOR_VALUE: {
68 const auto x = input(0);
69 const auto y = input(1);
72 output.access(0, 0) = (1.0 - y)/4.0;
73 output.access(0, 1) = 0.0;
75 output.access(1, 0) = 0.0;
76 output.access(1, 1) = (1.0 + x)/4.0;
78 output.access(2, 0) = -(1.0 + y)/4.0;
79 output.access(2, 1) = 0.0;
81 output.access(3, 0) = 0.0;
82 output.access(3, 1) = -(1.0 - x)/4.0;
87 output.access(0) = 0.25;
88 output.access(1) = 0.25;
89 output.access(2) = 0.25;
90 output.access(3) = 0.25;
94 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
95 opType != OPERATOR_CURL,
96 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::Serial::getValues) operator is not supported");
101 template<
typename SpT,
102 typename outputValueValueType,
class ...outputValueProperties,
103 typename inputPointValueType,
class ...inputPointProperties>
105 Basis_HCURL_QUAD_I1_FEM::
106 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
107 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
108 const EOperator operatorType ) {
109 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
110 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
111 typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
114 const auto loopSize = inputPoints.extent(0);
115 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
117 switch (operatorType) {
118 case OPERATOR_VALUE: {
119 typedef Functor<outputValueViewType, inputPointViewType, OPERATOR_VALUE> FunctorType;
120 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
123 case OPERATOR_CURL: {
124 typedef Functor<outputValueViewType, inputPointViewType, OPERATOR_CURL> FunctorType;
125 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
129 INTREPID2_TEST_FOR_EXCEPTION( operatorType == OPERATOR_DIV, std::invalid_argument,
130 ">>> ERROR (Basis_HCURL_QUAD_I1_FEM::getValues): DIV is invalid operator for HCURL Basis Functions");
133 case OPERATOR_GRAD: {
134 INTREPID2_TEST_FOR_EXCEPTION( operatorType == OPERATOR_GRAD, std::invalid_argument,
135 ">>> ERROR (Basis_HCURL_QUAD_I1_FEM::getValues): GRAD is invalid operator for HCURL Basis Functions");
148 INTREPID2_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_D1) ||
149 (operatorType == OPERATOR_D2) ||
150 (operatorType == OPERATOR_D3) ||
151 (operatorType == OPERATOR_D4) ||
152 (operatorType == OPERATOR_D5) ||
153 (operatorType == OPERATOR_D6) ||
154 (operatorType == OPERATOR_D7) ||
155 (operatorType == OPERATOR_D8) ||
156 (operatorType == OPERATOR_D9) ||
157 (operatorType == OPERATOR_D10),
158 std::invalid_argument,
159 ">>> ERROR (Basis_HCURL_QUAD_I1_FEM::getValues): Invalid operator type");
163 INTREPID2_TEST_FOR_EXCEPTION( (operatorType != OPERATOR_VALUE) &&
164 (operatorType != OPERATOR_GRAD) &&
165 (operatorType != OPERATOR_CURL) &&
166 (operatorType != OPERATOR_DIV) &&
167 (operatorType != OPERATOR_D1) &&
168 (operatorType != OPERATOR_D2) &&
169 (operatorType != OPERATOR_D3) &&
170 (operatorType != OPERATOR_D4) &&
171 (operatorType != OPERATOR_D5) &&
172 (operatorType != OPERATOR_D6) &&
173 (operatorType != OPERATOR_D7) &&
174 (operatorType != OPERATOR_D8) &&
175 (operatorType != OPERATOR_D9) &&
176 (operatorType != OPERATOR_D10),
177 std::invalid_argument,
178 ">>> ERROR (Basis_HCURL_QUAD_I1_FEM::getValues): Invalid operator type");
187 template<
typename SpT,
typename OT,
typename PT>
190 this->basisCardinality_ = 4;
191 this->basisDegree_ = 1;
192 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
193 this->basisType_ = BASIS_FEM_DEFAULT;
194 this->basisCoordinates_ = COORDINATES_CARTESIAN;
195 this->functionSpace_ = FUNCTION_SPACE_HCURL;
200 const ordinal_type tagSize = 4;
201 const ordinal_type posScDim = 0;
202 const ordinal_type posScOrd = 1;
203 const ordinal_type posDfOrd = 2;
206 ordinal_type tags[16] = { 1, 0, 0, 1,
216 this->setOrdinalTagData(this->tagToOrdinal_,
219 this->basisCardinality_,
227 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
228 dofCoords(
"dofCoordsHost", this->basisCardinality_,this->basisCellTopology_.getDimension());
230 dofCoords(0,0) = 0.0; dofCoords(0,1) = -1.0;
231 dofCoords(1,0) = 1.0; dofCoords(1,1) = 0.0;
232 dofCoords(2,0) = 0.0; dofCoords(2,1) = 1.0;
233 dofCoords(3,0) = -1.0; dofCoords(3,1) = 0.0;
235 this->dofCoords_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoords);
236 Kokkos::deep_copy(this->dofCoords_, dofCoords);
240 Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
241 dofCoeffs(
"dofCoeffsHost", this->basisCardinality_,this->basisCellTopology_.getDimension());
244 dofCoeffs(0,0) = 2.0; dofCoeffs(0,1) = 0.0;
245 dofCoeffs(1,0) = 0.0; dofCoeffs(1,1) = 2.0;
246 dofCoeffs(2,0) = -2.0; dofCoeffs(2,1) = 0.0;
247 dofCoeffs(3,0) = 0.0; dofCoeffs(3,1) = -2.0;
249 this->dofCoeffs_ = Kokkos::create_mirror_view(
typename SpT::memory_space(), dofCoeffs);
250 Kokkos::deep_copy(this->dofCoeffs_, dofCoeffs);
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Basis_HCURL_QUAD_I1_FEM()
Constructor.