Intrepid2
Intrepid2_HCURL_WEDGE_I1_FEMDef.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
49 #ifndef __INTREPID2_HCURL_WEDGE_I1_FEM_DEF_HPP__
50 #define __INTREPID2_HCURL_WEDGE_I1_FEM_DEF_HPP__
51 
52 namespace Intrepid2 {
53 
54  // -------------------------------------------------------------------------------------
55  namespace Impl {
56 
57  template<EOperator opType>
58  template<typename outputViewType,
59  typename inputViewType>
60  KOKKOS_INLINE_FUNCTION
61  void
62  Basis_HCURL_WEDGE_I1_FEM::Serial<opType>::
63  getValues( outputViewType output,
64  const inputViewType input ) {
65  switch (opType) {
66  case OPERATOR_VALUE: {
67  const auto x = input(0);
68  const auto y = input(1);
69  const auto z = input(2);
70 
71  // outputValues is a rank-3 array with dimensions (basisCardinality_, dim0, spaceDim)
72  output.access(0, 0) = (1.0 - z)*(1.0 - y)/2.0;
73  output.access(0, 1) = x*(1.0 - z)/2.0;
74  output.access(0, 2) = 0.0;
75 
76  output.access(1, 0) = y*(z - 1.0)/2.0;
77  output.access(1, 1) = x*(1.0 - z)/2.0;
78  output.access(1, 2) = 0.0;
79 
80  output.access(2, 0) = y*(z - 1.0)/2.0;
81  output.access(2, 1) = (1.0 - x)*(z - 1.0)/2.0;
82  output.access(2, 2) = 0.0;
83 
84  output.access(3, 0) = (1.0 - y)*(1.0 + z)/2.0;
85  output.access(3, 1) = x*(1.0 + z)/2.0;
86  output.access(3, 2) = 0.0;
87 
88  output.access(4, 0) =-y*(1.0 + z)/2.0;
89  output.access(4, 1) = x*(1.0 + z)/2.0;
90  output.access(4, 2) = 0.0;
91 
92  output.access(5, 0) = -y*(1.0 + z)/2.0;
93  output.access(5, 1) = (x - 1.0)*(1.0 + z)/2.0;
94  output.access(5, 2) = 0.0;
95 
96  output.access(6, 0) = 0.0;
97  output.access(6, 1) = 0.0;
98  output.access(6, 2) = (1.0 - x - y)/2.0;
99 
100  output.access(7, 0) = 0.0;
101  output.access(7, 1) = 0.0;
102  output.access(7, 2) = x/2.0;
103 
104  output.access(8, 0) = 0.0;
105  output.access(8, 1) = 0.0;
106  output.access(8, 2) = y/2.0;
107  break;
108  }
109  case OPERATOR_CURL: {
110  const auto x = input(0);
111  const auto y = input(1);
112  const auto z = input(2);
113  // outputValues is a rank-3 array with dimensions (basisCardinality_, dim0, spaceDim)
114  output.access(0, 0) = x/2.0;
115  output.access(0, 1) = (y - 1.0)/2.0;
116  output.access(0, 2) = 1.0 - z;
117 
118  output.access(1, 0) = x/2.0;
119  output.access(1, 1) = y/2.0;
120  output.access(1, 2) = 1.0 - z;
121 
122  output.access(2, 0) = (x - 1.0)/2.0;
123  output.access(2, 1) = y/2.0;
124  output.access(2, 2) = 1.0 - z;
125 
126  output.access(3, 0) = -x/2.0;
127  output.access(3, 1) = (1.0 - y)/2.0;
128  output.access(3, 2) = 1.0 + z;
129 
130  output.access(4, 0) = -x/2.0;
131  output.access(4, 1) = -y/2.0;
132  output.access(4, 2) = 1.0 + z;
133 
134  output.access(5, 0) = (1.0 - x)/2.0;
135  output.access(5, 1) = -y/2.0;
136  output.access(5, 2) = 1.0 + z;
137 
138  output.access(6, 0) =-0.5;
139  output.access(6, 1) = 0.5;
140  output.access(6, 2) = 0.0;
141 
142  output.access(7, 0) = 0.0;
143  output.access(7, 1) =-0.5;
144  output.access(7, 2) = 0.0;
145 
146  output.access(8, 0) = 0.5;
147  output.access(8, 1) = 0.0;
148  output.access(8, 2) = 0.0;
149  break;
150  }
151  default: {
152  INTREPID2_TEST_FOR_ABORT( opType != OPERATOR_VALUE &&
153  opType != OPERATOR_CURL,
154  ">>> ERROR: (Intrepid2::Basis_HCURL_WEDGE_I1_FEM::Serial::getValues) operator is not supported");
155 
156  }
157  }
158  }
159 
160  template<typename SpT,
161  typename outputValueValueType, class ...outputValueProperties,
162  typename inputPointValueType, class ...inputPointProperties>
163  void Basis_HCURL_WEDGE_I1_FEM::
164  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
165  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
166  const EOperator operatorType ) {
167  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
168  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
169  typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
170 
171  // Number of evaluation points = dim 0 of inputPoints
172  const auto loopSize = inputPoints.extent(0);
173  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
174 
175  switch (operatorType) {
176  case OPERATOR_VALUE: {
177  typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_VALUE> FunctorType;
178  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
179  break;
180  }
181  case OPERATOR_CURL: {
182  typedef Functor<outputValueViewType,inputPointViewType,OPERATOR_CURL> FunctorType;
183  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints) );
184  break;
185  }
186  case OPERATOR_DIV: {
187  INTREPID2_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_DIV), std::invalid_argument,
188  ">>> ERROR (Basis_HCURL_WEDGE_I1_FEM): DIV is invalid operator for HCURL Basis Functions");
189  break;
190  }
191  case OPERATOR_GRAD: {
192  INTREPID2_TEST_FOR_EXCEPTION( (operatorType == OPERATOR_GRAD), std::invalid_argument,
193  ">>> ERROR (Basis_HCURL_WEDGE_I1_FEM): GRAD is invalid operator for HCURL Basis Functions");
194  break;
195  }
196  case OPERATOR_D1:
197  case OPERATOR_D2:
198  case OPERATOR_D3:
199  case OPERATOR_D4:
200  case OPERATOR_D5:
201  case OPERATOR_D6:
202  case OPERATOR_D7:
203  case OPERATOR_D8:
204  case OPERATOR_D9:
205  case OPERATOR_D10: {
206  INTREPID2_TEST_FOR_EXCEPTION( ( (operatorType == OPERATOR_D1) ||
207  (operatorType == OPERATOR_D2) ||
208  (operatorType == OPERATOR_D3) ||
209  (operatorType == OPERATOR_D4) ||
210  (operatorType == OPERATOR_D5) ||
211  (operatorType == OPERATOR_D6) ||
212  (operatorType == OPERATOR_D7) ||
213  (operatorType == OPERATOR_D8) ||
214  (operatorType == OPERATOR_D9) ||
215  (operatorType == OPERATOR_D10) ),
216  std::invalid_argument,
217  ">>> ERROR (Basis_HCURL_WEDGE_I1_FEM): Invalid operator type");
218  break;
219  }
220  default: {
221  INTREPID2_TEST_FOR_EXCEPTION( !(Intrepid2::isValidOperator(operatorType)),
222  std::invalid_argument,
223  ">>> ERROR (Basis_HCURL_WEDGE_I1_FEM): Invalid operator type");
224  }
225  }
226  }
227 
228  }
229 
230  template<typename SpT, typename OT, typename PT>
233  this->basisCardinality_ = 9;
234  this->basisDegree_ = 1;
235  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Wedge<6> >() );
236  this->basisType_ = BASIS_FEM_DEFAULT;
237  this->basisCoordinates_ = COORDINATES_CARTESIAN;
238 
239 
240  // initialize tags
241  {
242  // Basis-dependent intializations
243  const ordinal_type tagSize = 4; // size of DoF tag
244  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
245  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
246  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
247 
248  // An array with local DoF tags assigned to basis functions, in the order of their local enumeration
249  ordinal_type tags[36] = {
250  1, 0, 0, 1,
251  1, 1, 0, 1,
252  1, 2, 0, 1,
253  1, 3, 0, 1,
254  1, 4, 0, 1,
255  1, 5, 0, 1,
256  1, 6, 0, 1,
257  1, 7, 0, 1,
258  1, 8, 0, 1 };
259 
260  //host tags
261  ordinal_type_array_1d_host tagView(&tags[0],36);
262 
263  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
264  this->setOrdinalTagData(this->tagToOrdinal_,
265  this->ordinalToTag_,
266  tagView,
267  this->basisCardinality_,
268  tagSize,
269  posScDim,
270  posScOrd,
271  posDfOrd);
272  }
273  // dofCoords on host and create its mirror view to device
274  Kokkos::DynRankView<typename scalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
275  dofCoords("dofCoordsHost", this->basisCardinality_,this->basisCellTopology_.getDimension());
276 
277  dofCoords(0,0) = 0.5; dofCoords(0,1) = 0.0; dofCoords(0,2) = -1.0;
278  dofCoords(1,0) = 0.5; dofCoords(1,1) = 0.5; dofCoords(1,2) = -1.0;
279  dofCoords(2,0) = 0.0; dofCoords(2,1) = 0.5; dofCoords(2,2) = -1.0;
280  dofCoords(3,0) = 0.5; dofCoords(3,1) = 0.0; dofCoords(3,2) = 1.0;
281  dofCoords(4,0) = 0.5; dofCoords(4,1) = 0.5; dofCoords(4,2) = 1.0;
282  dofCoords(5,0) = 0.0; dofCoords(5,1) = 0.5; dofCoords(5,2) = 1.0;
283  dofCoords(6,0) = 0.0; dofCoords(6,1) = 0.0; dofCoords(6,2) = 0.0;
284  dofCoords(7,0) = 1.0; dofCoords(7,1) = 0.0; dofCoords(7,2) = 0.0;
285  dofCoords(8,0) = 0.0; dofCoords(8,1) = 1.0; dofCoords(8,2) = 0.0;
286 
287  this->dofCoords_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoords);
288  Kokkos::deep_copy(this->dofCoords_, dofCoords);
289 
290  }
291 
292 }// namespace Intrepid2
293 #endif
Kokkos::View< ordinal_type *,typename ExecSpaceType::array_layout, Kokkos::HostSpace > ordinal_type_array_1d_host
View type for 1d host array.