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