49#ifndef __INTREPID2_HGRAD_QUAD_C1_FEM_DEF_HPP__
50#define __INTREPID2_HGRAD_QUAD_C1_FEM_DEF_HPP__
58 template<EOperator opType>
59 template<
typename OutputViewType,
60 typename inputViewType>
61 KOKKOS_INLINE_FUNCTION
65 const inputViewType input ) {
67 case OPERATOR_VALUE : {
68 const auto x = input(0);
69 const auto y = input(1);
72 output.access(0) = (1.0 - x)*(1.0 - y)/4.0;
73 output.access(1) = (1.0 + x)*(1.0 - y)/4.0;
74 output.access(2) = (1.0 + x)*(1.0 + y)/4.0;
75 output.access(3) = (1.0 - x)*(1.0 + y)/4.0;
78 case OPERATOR_GRAD : {
79 const auto x = input(0);
80 const auto y = input(1);
83 output.access(0, 0) = -(1.0 - y)/4.0;
84 output.access(0, 1) = -(1.0 - x)/4.0;
86 output.access(1, 0) = (1.0 - y)/4.0;
87 output.access(1, 1) = -(1.0 + x)/4.0;
89 output.access(2, 0) = (1.0 + y)/4.0;
90 output.access(2, 1) = (1.0 + x)/4.0;
92 output.access(3, 0) = -(1.0 + y)/4.0;
93 output.access(3, 1) = (1.0 - x)/4.0;
96 case OPERATOR_CURL : {
97 const auto x = input(0);
98 const auto y = input(1);
101 output.access(0, 0) = -(1.0 - x)/4.0;
102 output.access(0, 1) = (1.0 - y)/4.0;
104 output.access(1, 0) = -(1.0 + x)/4.0;
105 output.access(1, 1) = -(1.0 - y)/4.0;
107 output.access(2, 0) = (1.0 + x)/4.0;
108 output.access(2, 1) = -(1.0 + y)/4.0;
110 output.access(3, 0) = (1.0 - x)/4.0;
111 output.access(3, 1) = (1.0 + y)/4.0;
116 output.access(0, 0) = 0.0;
117 output.access(0, 1) = 0.25;
118 output.access(0, 2) = 0.0;
120 output.access(1, 0) = 0.0;
121 output.access(1, 1) = -0.25;
122 output.access(1, 2) = 0.0;
124 output.access(2, 0) = 0.0;
125 output.access(2, 1) = 0.25;
126 output.access(2, 2) = 0.0;
128 output.access(3, 0) = 0.0;
129 output.access(3, 1) = -0.25;
130 output.access(3, 2) = 0.0;
133 case OPERATOR_MAX : {
134 const ordinal_type jend = output.extent(1);
135 const ordinal_type iend = output.extent(0);
137 for (ordinal_type j=0;j<jend;++j)
138 for (ordinal_type i=0;i<iend;++i)
139 output.access(i, j) = 0.0;
143 INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
144 opType != OPERATOR_GRAD &&
145 opType != OPERATOR_CURL &&
146 opType != OPERATOR_D2 &&
147 opType != OPERATOR_MAX,
148 ">>> ERROR: (Intrepid2::Basis_HGRAD_QUAD_C1_FEM::Serial::getValues) operator is not supported");
154 template<
typename DT,
155 typename outputValueValueType,
class ...outputValueProperties,
156 typename inputPointValueType,
class ...inputPointProperties>
158 Basis_HGRAD_QUAD_C1_FEM::
159 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
160 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
161 const EOperator operatorType ) {
162 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
163 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
164 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
167 const auto loopSize = inputPoints.extent(0);
168 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
170 switch (operatorType) {
172 case OPERATOR_VALUE: {
174 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
181 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
184 case OPERATOR_CURL: {
186 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
190 INTREPID2_TEST_FOR_EXCEPTION( operatorType == OPERATOR_DIV, std::invalid_argument,
191 ">>> ERROR (Basis_HGRAD_QUAD_C1_FEM): DIV is invalid operator for rank-0 (scalar) functions in 2D");
196 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
208 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
212 INTREPID2_TEST_FOR_EXCEPTION( !Intrepid2::isValidOperator(operatorType), std::invalid_argument,
213 ">>> ERROR (Basis_HGRAD_QUAD_C1_FEM): Invalid operator type");
220 template<
typename DT,
typename OT,
typename PT>
225 this->
basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
233 const ordinal_type tagSize = 4;
234 const ordinal_type posScDim = 0;
235 const ordinal_type posScOrd = 1;
236 const ordinal_type posDfOrd = 2;
239 ordinal_type tags[16] = { 0, 0, 0, 1,
245 OrdinalTypeArray1DHost tagView(&tags[0], 16);
267 Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
270 dofCoords(0,0) = -1.0; dofCoords(0,1) = -1.0;
271 dofCoords(1,0) = 1.0; dofCoords(1,1) = -1.0;
272 dofCoords(2,0) = 1.0; dofCoords(2,1) = 1.0;
273 dofCoords(3,0) = -1.0; dofCoords(3,1) = 1.0;
275 this->
dofCoords_ = Kokkos::create_mirror_view(
typename DT::memory_space(), dofCoords);
276 Kokkos::deep_copy(this->
dofCoords_, dofCoords);
Basis_HGRAD_QUAD_C1_FEM()
Constructor.
ECoordinates basisCoordinates_
ordinal_type basisDegree_
void setOrdinalTagData(OrdinalTypeView3D &tagToOrdinal, OrdinalTypeView2D &ordinalToTag, const OrdinalTypeView1D tags, const ordinal_type basisCard, const ordinal_type tagSize, const ordinal_type posScDim, const ordinal_type posScOrd, const ordinal_type posDfOrd)
OrdinalTypeArray2DHost ordinalToTag_
Kokkos::DynRankView< scalarType, DeviceType > dofCoords_
ordinal_type basisCardinality_
OrdinalTypeArray3DHost tagToOrdinal_
shards::CellTopology basisCellTopology_
EFunctionSpace functionSpace_
See Intrepid2::Basis_HGRAD_QUAD_C1_FEM.
See Intrepid2::Basis_HGRAD_QUAD_C1_FEM.