Intrepid2
Intrepid2_HGRAD_LINE_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 
49 #ifndef __INTREPID2_HGRAD_LINE_CN_FEM_DEF_HPP__
50 #define __INTREPID2_HGRAD_LINE_CN_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  typename workViewType,
61  typename vinvViewType>
62  KOKKOS_INLINE_FUNCTION
63  void
64  Basis_HGRAD_LINE_Cn_FEM::Serial<opType>::
65  getValues( OutputViewType output,
66  const inputViewType input,
67  workViewType work,
68  const vinvViewType vinv,
69  const ordinal_type operatorDn ) {
70  ordinal_type opDn = operatorDn;
71 
72  const ordinal_type card = vinv.extent(0);
73  const ordinal_type npts = input.extent(0);
74 
75  const ordinal_type order = card - 1;
76  const double alpha = 0.0, beta = 0.0;
77 
78  typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
79  auto vcprop = Kokkos::common_view_alloc_prop(work);
80 
81  switch (opType) {
82  case OPERATOR_VALUE: {
83  viewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts);
84 
85  Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
86  Serial<opType>::getValues(phis, input, order, alpha, beta);
87 
88  for (ordinal_type i=0;i<card;++i)
89  for (ordinal_type j=0;j<npts;++j) {
90  output.access(i,j) = 0.0;
91  for (ordinal_type k=0;k<card;++k)
92  output.access(i,j) += vinv(k,i)*phis.access(k,j);
93  }
94  break;
95  }
96  case OPERATOR_GRAD:
97  case OPERATOR_D1:
98  case OPERATOR_D2:
99  case OPERATOR_D3:
100  case OPERATOR_D4:
101  case OPERATOR_D5:
102  case OPERATOR_D6:
103  case OPERATOR_D7:
104  case OPERATOR_D8:
105  case OPERATOR_D9:
106  case OPERATOR_D10:
107  opDn = getOperatorOrder(opType);
108  case OPERATOR_Dn: {
109  // dkcard is always 1 for 1D element
110  const ordinal_type dkcard = 1;
111  viewType phis(Kokkos::view_wrap(work.data(), vcprop), card, npts, dkcard);
112  Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
113  Serial<opType>::getValues(phis, input, order, alpha, beta, opDn);
114 
115  for (ordinal_type i=0;i<card;++i)
116  for (ordinal_type j=0;j<npts;++j)
117  for (ordinal_type k=0;k<dkcard;++k) {
118  output.access(i,j,k) = 0.0;
119  for (ordinal_type l=0;l<card;++l)
120  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
121  }
122  break;
123  }
124  default: {
125  INTREPID2_TEST_FOR_ABORT( true,
126  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM::Serial::getValues) operator is not supported." );
127  }
128  }
129  }
130 
131 
132  template<typename DT, ordinal_type numPtsPerEval,
133  typename outputValueValueType, class ...outputValueProperties,
134  typename inputPointValueType, class ...inputPointProperties,
135  typename vinvValueType, class ...vinvProperties>
136  void
137  Basis_HGRAD_LINE_Cn_FEM::
138  getValues( const typename DT::execution_space& space,
139  Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
140  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
141  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
142  const EOperator operatorType ) {
143  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
144  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
145  typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
146  typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
147 
148  // loopSize corresponds to cardinality
149  const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
150  const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
151  const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
152  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(space, 0, loopSize);
153 
154  typedef typename inputPointViewType::value_type inputPointType;
155 
156  const ordinal_type cardinality = outputValues.extent(0);
157 
158  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
159  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
160  workViewType work(Kokkos::view_alloc(space, "Basis_HGRAD_LINE_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
161 
162  switch (operatorType) {
163  case OPERATOR_VALUE: {
164  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
165  OPERATOR_VALUE,numPtsPerEval> FunctorType;
166  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
167  break;
168  }
169  case OPERATOR_GRAD:
170  case OPERATOR_D1:
171  case OPERATOR_D2:
172  case OPERATOR_D3:
173  case OPERATOR_D4:
174  case OPERATOR_D5:
175  case OPERATOR_D6:
176  case OPERATOR_D7:
177  case OPERATOR_D8:
178  case OPERATOR_D9:
179  case OPERATOR_D10: {
180  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,workViewType,
181  OPERATOR_Dn,numPtsPerEval> FunctorType;
182  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work,
183  getOperatorOrder(operatorType)) );
184  break;
185  }
186  default: {
187  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
188  ">>> ERROR (Basis_HGRAD_LINE_Cn_FEM): Operator type not implemented" );
189  //break; commented out because this always throws
190  }
191  }
192  }
193  }
194 
195  // -------------------------------------------------------------------------------------
196  template<typename DT, typename OT, typename PT>
198  Basis_HGRAD_LINE_Cn_FEM( const ordinal_type order,
199  const EPointType pointType ) {
200  this->pointType_ = pointType;
201  this->basisCardinality_ = order+1;
202  this->basisDegree_ = order;
203  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Line<2> >() );
204  this->basisType_ = BASIS_FEM_LAGRANGIAN;
205  this->basisCoordinates_ = COORDINATES_CARTESIAN;
206  this->functionSpace_ = FUNCTION_SPACE_HGRAD;
207 
208  const ordinal_type card = this->basisCardinality_;
209 
210  // points are computed in the host and will be copied
211  Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
212  dofCoords("Hgrad::Line::Cn::dofCoords", card, 1);
213 
214  //Default is Equispaced
215  auto pointT = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
216 
217  switch (pointT) {
218  case POINTTYPE_EQUISPACED:
219  case POINTTYPE_WARPBLEND: {
220  // lattice ordering
221  {
222  const ordinal_type offset = 0;
223  PointTools::getLattice( dofCoords,
224  this->basisCellTopology_,
225  order, offset,
226  pointT );
227 
228  }
229  break;
230  }
231  default: {
232  INTREPID2_TEST_FOR_EXCEPTION( !isValidPointType(pointT),
233  std::invalid_argument ,
234  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM) invalid pointType." );
235  }
236  }
237 
238  this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
239  Kokkos::deep_copy(this->dofCoords_, dofCoords);
240 
241  // form Vandermonde matrix; actually, this is the transpose of the VDM,
242  // this matrix is used in LAPACK so it should be column major and left layout
243  const ordinal_type lwork = card*card;
244  Kokkos::DynRankView<typename ScalarViewType::value_type,Kokkos::LayoutLeft,Kokkos::HostSpace>
245  vmat("Hgrad::Line::Cn::vmat", card, card),
246  work("Hgrad::Line::Cn::work", lwork),
247  ipiv("Hgrad::Line::Cn::ipiv", card);
248 
249  const double alpha = 0.0, beta = 0.0;
250  Impl::Basis_HGRAD_LINE_Cn_FEM_JACOBI::
251  getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>
252  (typename Kokkos::HostSpace::execution_space{}, vmat, dofCoords, order, alpha, beta, OPERATOR_VALUE);
253 
254  ordinal_type info = 0;
255  Teuchos::LAPACK<ordinal_type,typename ScalarViewType::value_type> lapack;
256 
257  lapack.GETRF(card, card,
258  vmat.data(), vmat.stride_1(),
259  (ordinal_type*)ipiv.data(),
260  &info);
261 
262  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
263  std::runtime_error ,
264  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM) lapack.GETRF returns nonzero info." );
265 
266  lapack.GETRI(card,
267  vmat.data(), vmat.stride_1(),
268  (ordinal_type*)ipiv.data(),
269  work.data(), lwork,
270  &info);
271 
272  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
273  std::runtime_error ,
274  ">>> ERROR: (Intrepid2::Basis_HGRAD_LINE_Cn_FEM) lapack.GETRI returns nonzero info." );
275 
276  // create host mirror
277  Kokkos::DynRankView<typename ScalarViewType::value_type,typename DT::execution_space::array_layout,Kokkos::HostSpace>
278  vinv("Hgrad::Line::Cn::vinv", card, card);
279 
280  for (ordinal_type i=0;i<card;++i)
281  for (ordinal_type j=0;j<card;++j)
282  vinv(i,j) = vmat(j,i);
283 
284  this->vinv_ = Kokkos::create_mirror_view(typename DT::memory_space(), vinv);
285  Kokkos::deep_copy(this->vinv_ , vinv);
286 
287  // initialize tags
288  {
289  // Basis-dependent initializations
290  const ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
291  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
292  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
293  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
294 
295  // Note: the only reason why equispaced can't support higher order than Parameters::MaxOrder appears to be the fact that the tags below get stored into a fixed-length array.
296  // TODO: relax the maximum order requirement by setting up tags in a different container, perhaps directly into an OrdinalTypeArray1DHost (tagView, below). (As of this writing (1/25/22), looks like other nodal bases do this in a similar way -- those should be fixed at the same time; maybe search for Parameters::MaxOrder.)
297  INTREPID2_TEST_FOR_EXCEPTION( order > Parameters::MaxOrder, std::invalid_argument, "polynomial order exceeds the max supported by this class");
298  ordinal_type tags[Parameters::MaxOrder+1][4];
299 
300  // lattice order
301  {
302  const auto v0 = 0;
303  tags[v0][0] = 0; // vertex dof
304  tags[v0][1] = 0; // vertex id
305  tags[v0][2] = 0; // local dof id
306  tags[v0][3] = 1; // total number of dofs in this vertex
307 
308  const ordinal_type iend = card - 2;
309  for (ordinal_type i=0;i<iend;++i) {
310  const auto e = i + 1;
311  tags[e][0] = 1; // edge dof
312  tags[e][1] = 0; // edge id
313  tags[e][2] = i; // local dof id
314  tags[e][3] = iend; // total number of dofs in this edge
315  }
316 
317  const auto v1 = card -1;
318  tags[v1][0] = 0; // vertex dof
319  tags[v1][1] = 1; // vertex id
320  tags[v1][2] = 0; // local dof id
321  tags[v1][3] = 1; // total number of dofs in this vertex
322  }
323 
324  // topological order
325  // {
326  // tags[0][0] = 0; // vertex dof
327  // tags[0][1] = 0; // vertex id
328  // tags[0][2] = 0; // local dof id
329  // tags[0][3] = 1; // total number of dofs in this vertex
330 
331  // tags[1][0] = 0; // vertex dof
332  // tags[1][1] = 1; // vertex id
333  // tags[1][2] = 0; // local dof id
334  // tags[1][3] = 1; // total number of dofs in this vertex
335 
336  // const ordinal_type iend = card - 2;
337  // for (ordinal_type i=0;i<iend;++i) {
338  // const auto ii = i + 2;
339  // tags[ii][0] = 1; // edge dof
340  // tags[ii][1] = 0; // edge id
341  // tags[ii][2] = i; // local dof id
342  // tags[ii][3] = iend; // total number of dofs in this edge
343  // }
344  // }
345 
346 
347  OrdinalTypeArray1DHost tagView(&tags[0][0], card*4);
348 
349  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
350  // tags are constructed on host
351  this->setOrdinalTagData(this->tagToOrdinal_,
352  this->ordinalToTag_,
353  tagView,
354  this->basisCardinality_,
355  tagSize,
356  posScDim,
357  posScOrd,
358  posDfOrd);
359  }
360  }
361 
362 }// namespace Intrepid2
363 
364 #endif
365 
366 
367 
368 
369 
370 
371 
372 
373 
374 
375 
376 
377 
378 
379 
static void getLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference simplex, quadrilateral or hexahedron (cu...
Basis_HGRAD_LINE_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
Kokkos::View< ordinal_type *, typename ExecutionSpace::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.