Intrepid2
Intrepid2_HVOL_QUAD_Cn_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 
48 #ifndef __INTREPID2_HVOL_QUAD_CN_FEM_DEF_HPP__
49 #define __INTREPID2_HVOL_QUAD_CN_FEM_DEF_HPP__
50 
51 namespace Intrepid2 {
52 
53  // -------------------------------------------------------------------------------------
54  namespace Impl {
55 
56  template<EOperator opType>
57  template<typename OutputViewType,
58  typename inputViewType,
59  typename workViewType,
60  typename vinvViewType>
61  KOKKOS_INLINE_FUNCTION
62  void
63  Basis_HVOL_QUAD_Cn_FEM::Serial<opType>::
64  getValues( OutputViewType output,
65  const inputViewType input,
66  workViewType work,
67  const vinvViewType vinv,
68  const ordinal_type operatorDn ) {
69  ordinal_type opDn = operatorDn;
70 
71  const ordinal_type cardLine = vinv.extent(0);
72  const ordinal_type npts = input.extent(0);
73 
74  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
75  const auto input_x = Kokkos::subview(input, Kokkos::ALL(), range_type(0,1));
76  const auto input_y = Kokkos::subview(input, Kokkos::ALL(), range_type(1,2));
77 
78  const int dim_s = get_dimension_scalar(work);
79  auto ptr0 = work.data();
80  auto ptr1 = work.data()+cardLine*npts*dim_s;
81  auto ptr2 = work.data()+2*cardLine*npts*dim_s;
82 
83  typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
84  auto vcprop = Kokkos::common_view_alloc_prop(work);
85 
86  switch (opType) {
87  case OPERATOR_VALUE: {
88  viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
89  viewType output_x(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
90  viewType output_y(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
91 
92  Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
93  getValues(output_x, input_x, work_line, vinv);
94 
95  Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
96  getValues(output_y, input_y, work_line, vinv);
97 
98  // tensor product
99  ordinal_type idx = 0;
100  for (ordinal_type j=0;j<cardLine;++j) // y
101  for (ordinal_type i=0;i<cardLine;++i,++idx) // x
102  for (ordinal_type k=0;k<npts;++k)
103  output.access(idx,k) = output_x.access(i,k)*output_y.access(j,k);
104  break;
105  }
106  case OPERATOR_GRAD:
107  case OPERATOR_D1:
108  case OPERATOR_D2:
109  case OPERATOR_D3:
110  case OPERATOR_D4:
111  case OPERATOR_D5:
112  case OPERATOR_D6:
113  case OPERATOR_D7:
114  case OPERATOR_D8:
115  case OPERATOR_D9:
116  case OPERATOR_D10:
117  opDn = getOperatorOrder(opType);
118  case OPERATOR_Dn: {
119  const auto dkcard = opDn + 1;
120  for (auto l=0;l<dkcard;++l) {
121  viewType work_line(Kokkos::view_wrap(ptr0, vcprop), cardLine, npts);
122 
123  viewType output_x, output_y;
124 
125  const auto mult_x = opDn - l;
126  const auto mult_y = l;
127 
128  if (mult_x) {
129  output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts, 1);
130  Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
131  getValues(output_x, input_x, work_line, vinv, mult_x);
132  } else {
133  output_x = viewType(Kokkos::view_wrap(ptr1, vcprop), cardLine, npts);
134  Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
135  getValues(output_x, input_x, work_line, vinv);
136  }
137 
138  if (mult_y) {
139  output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts, 1);
140  Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_Dn>::
141  getValues(output_y, input_y, work_line, vinv, mult_y);
142  } else {
143  output_y = viewType(Kokkos::view_wrap(ptr2, vcprop), cardLine, npts);
144  Impl::Basis_HVOL_LINE_Cn_FEM::Serial<OPERATOR_VALUE>::
145  getValues(output_y, input_y, work_line, vinv);
146  }
147 
148  // tensor product (extra dimension of ouput x and y are ignored)
149  ordinal_type idx = 0;
150  for (ordinal_type j=0;j<cardLine;++j) // y
151  for (ordinal_type i=0;i<cardLine;++i,++idx) // x
152  for (ordinal_type k=0;k<npts;++k)
153  output.access(idx,k,l) = output_x.access(i,k,0)*output_y.access(j,k,0);
154  }
155  break;
156  }
157  default: {
158  INTREPID2_TEST_FOR_ABORT( true,
159  ">>> ERROR: (Intrepid2::Basis_HVOL_QUAD_Cn_FEM::Serial::getValues) operator is not supported" );
160  }
161  }
162  }
163 
164  template<typename SpT, ordinal_type numPtsPerEval,
165  typename outputValueValueType, class ...outputValueProperties,
166  typename inputPointValueType, class ...inputPointProperties,
167  typename vinvValueType, class ...vinvProperties>
168  void
169  Basis_HVOL_QUAD_Cn_FEM::
170  getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
171  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
172  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
173  const EOperator operatorType ) {
174  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
175  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
176  typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
177  typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
178 
179  // loopSize corresponds to cardinality
180  const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
181  const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
182  const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
183  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
184 
185  typedef typename inputPointViewType::value_type inputPointType;
186 
187  const ordinal_type cardinality = outputValues.extent(0);
188  const ordinal_type cardLine = std::sqrt(cardinality);
189  const ordinal_type workSize = 3*cardLine;
190 
191  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
192  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
193  workViewType work(Kokkos::view_alloc("Basis_HVOL_QUAD_Cn_FEM::getValues::work", vcprop), workSize, inputPoints.extent(0));
194 
195  switch (operatorType) {
196  case OPERATOR_VALUE: {
197  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
198  OPERATOR_VALUE,numPtsPerEval> FunctorType;
199  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
200  break;
201  }
202  case OPERATOR_GRAD:
203  case OPERATOR_D1:
204  case OPERATOR_D2:
205  case OPERATOR_D3:
206  case OPERATOR_D4:
207  case OPERATOR_D5:
208  case OPERATOR_D6:
209  case OPERATOR_D7:
210  case OPERATOR_D8:
211  case OPERATOR_D9:
212  case OPERATOR_D10: {
213  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
214  OPERATOR_Dn,numPtsPerEval> FunctorType;
215  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
216  getOperatorOrder(operatorType)) );
217  break;
218  }
219  default: {
220  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
221  ">>> ERROR (Basis_HVOL_QUAD_Cn_FEM): Operator type not implemented" );
222  // break;commented out because exception
223  }
224  }
225  }
226  }
227 
228  // -------------------------------------------------------------------------------------
229  template<typename SpT, typename OT, typename PT>
231  Basis_HVOL_QUAD_Cn_FEM( const ordinal_type order,
232  const EPointType pointType ) {
233  // INTREPID2_TEST_FOR_EXCEPTION( !(pointType == POINTTYPE_EQUISPACED ||
234  // pointType == POINTTYPE_WARPBLEND), std::invalid_argument,
235  // ">>> ERROR (Basis_HVOL_QUAD_Cn_FEM): pointType must be either equispaced or warpblend." );
236 
237  // this should be in host
238  Basis_HVOL_LINE_Cn_FEM<SpT,OT,PT> lineBasis( order, pointType );
239  const auto cardLine = lineBasis.getCardinality();
240 
241  this->vinv_ = Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>("HVOL::Quad::Cn::vinv", cardLine, cardLine);
242  lineBasis.getVandermondeInverse(this->vinv_);
243 
244  this->basisCardinality_ = cardLine*cardLine;
245  this->basisDegree_ = order;
246  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Quadrilateral<4> >() );
247  this->basisType_ = BASIS_FEM_FIAT;
248  this->basisCoordinates_ = COORDINATES_CARTESIAN;
249  this->functionSpace_ = FUNCTION_SPACE_HVOL;
250 
251  // initialize tags
252  {
253  // Basis-dependent initializations
254  const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
255  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
256  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
257  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
258 
259  // An array with local DoF tags assigned to the basis functions, in the order of their local enumeration
260  constexpr ordinal_type maxCardLine = Parameters::MaxOrder + 1;
261  ordinal_type tags[maxCardLine*maxCardLine][4];
262 
263  {
264  ordinal_type idx = 0;
265  for (ordinal_type j=0;j<cardLine;++j) { // y
266  const auto tag_y = lineBasis.getDofTag(j);
267  for (ordinal_type i=0;i<cardLine;++i,++idx) { // x
268  const auto tag_x = lineBasis.getDofTag(i);
269 
270  // interior
271  tags[idx][0] = 2; // interior dof
272  tags[idx][1] = 0;
273  tags[idx][2] = tag_x(2) + tag_x(3)*tag_y(2); // local dof id
274  tags[idx][3] = tag_x(3)*tag_y(3); // total number of dofs in this vertex
275  }
276  }
277  }
278 
279  OrdinalTypeArray1DHost tagView(&tags[0][0], this->basisCardinality_*4);
280 
281  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
282  // tags are constructed on host
283  this->setOrdinalTagData(this->tagToOrdinal_,
284  this->ordinalToTag_,
285  tagView,
286  this->basisCardinality_,
287  tagSize,
288  posScDim,
289  posScOrd,
290  posDfOrd);
291  }
292 
293  // dofCoords on host and create its mirror view to device
294  Kokkos::DynRankView<typename ScalarViewType::value_type,typename SpT::array_layout,Kokkos::HostSpace>
295  dofCoordsHost("dofCoordsHost", this->basisCardinality_, this->basisCellTopology_.getDimension());
296 
297  Kokkos::DynRankView<typename ScalarViewType::value_type,SpT>
298  dofCoordsLine("dofCoordsLine", cardLine, 1);
299 
300  lineBasis.getDofCoords(dofCoordsLine);
301  auto dofCoordsLineHost = Kokkos::create_mirror_view(dofCoordsLine);
302  Kokkos::deep_copy(dofCoordsLineHost, dofCoordsLine);
303  {
304  ordinal_type idx = 0;
305  for (ordinal_type j=0;j<cardLine;++j) { // y
306  for (ordinal_type i=0;i<cardLine;++i,++idx) { // x
307  dofCoordsHost(idx,0) = dofCoordsLineHost(i,0);
308  dofCoordsHost(idx,1) = dofCoordsLineHost(j,0);
309  }
310  }
311  }
312 
313  this->dofCoords_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoordsHost);
314  Kokkos::deep_copy(this->dofCoords_, dofCoordsHost);
315  }
316 
317 }
318 
319 #endif
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Basis_HVOL_QUAD_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
ordinal_type getCardinality() const
Returns cardinality of the basis.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
const OrdinalTypeArrayStride1DHost getDofTag(const ordinal_type dofOrd) const
DoF ordinal to DoF tag lookup.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.