Intrepid2
Intrepid2_HVOL_TRI_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_TRI_CN_FEM_DEF_HPP__
49 #define __INTREPID2_HVOL_TRI_CN_FEM_DEF_HPP__
50 
52 
53 namespace Intrepid2 {
54 
55 // -------------------------------------------------------------------------------------
56 namespace Impl {
57 
58 template<EOperator opType>
59 template<typename OutputViewType,
60 typename inputViewType,
61 typename workViewType,
62 typename vinvViewType>
63 KOKKOS_INLINE_FUNCTION
64 void
65 Basis_HVOL_TRI_Cn_FEM::Serial<opType>::
66 getValues( OutputViewType output,
67  const inputViewType input,
68  workViewType work,
69  const vinvViewType vinv ) {
70 
71  constexpr ordinal_type spaceDim = 2;
72  const ordinal_type
73  card = vinv.extent(0),
74  npts = input.extent(0);
75 
76  // compute order
77  ordinal_type order = 0;
78  for (ordinal_type p=0;p<=Parameters::MaxOrder;++p) {
79  if (card == Intrepid2::getPnCardinality<spaceDim>(p) ) {
80  order = p;
81  break;
82  }
83  }
84 
85  typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
86  auto vcprop = Kokkos::common_view_alloc_prop(work);
87  auto ptr = work.data();
88 
89  switch (opType) {
90  case OPERATOR_VALUE: {
91  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
92  workViewType dummyView;
93 
94  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
95  Serial<opType>::getValues(phis, input, dummyView, order);
96 
97  for (ordinal_type i=0;i<card;++i)
98  for (ordinal_type j=0;j<npts;++j) {
99  output.access(i,j) = 0.0;
100  for (ordinal_type k=0;k<card;++k)
101  output.access(i,j) += vinv(k,i)*phis.access(k,j);
102  }
103  break;
104  }
105  case OPERATOR_GRAD:
106  case OPERATOR_D1: {
107  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
108  ptr += card*npts*spaceDim*get_dimension_scalar(work);
109  const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
110  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
111  Serial<opType>::getValues(phis, input, workView, order);
112 
113  for (ordinal_type i=0;i<card;++i)
114  for (ordinal_type j=0;j<npts;++j)
115  for (ordinal_type k=0;k<spaceDim;++k) {
116  output.access(i,j,k) = 0.0;
117  for (ordinal_type l=0;l<card;++l)
118  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
119  }
120  break;
121  }
122  case OPERATOR_D2:
123  case OPERATOR_D3:
124  case OPERATOR_D4:
125  case OPERATOR_D5:
126  case OPERATOR_D6:
127  case OPERATOR_D7:
128  case OPERATOR_D8:
129  case OPERATOR_D9:
130  case OPERATOR_D10: {
131  const ordinal_type dkcard = getDkCardinality<opType,spaceDim>(); //(orDn + 1);
132  const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
133  workViewType dummyView;
134 
135  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
136  Serial<opType>::getValues(phis, input, dummyView, order);
137 
138  for (ordinal_type i=0;i<card;++i)
139  for (ordinal_type j=0;j<npts;++j)
140  for (ordinal_type k=0;k<dkcard;++k) {
141  output.access(i,j,k) = 0.0;
142  for (ordinal_type l=0;l<card;++l)
143  output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
144  }
145  break;
146  }
147  default: {
148  INTREPID2_TEST_FOR_ABORT( true,
149  ">>> ERROR (Basis_HVOL_TRI_Cn_FEM): Operator type not implemented");
150  }
151  }
152 }
153 
154 template<typename SpT, ordinal_type numPtsPerEval,
155 typename outputValueValueType, class ...outputValueProperties,
156 typename inputPointValueType, class ...inputPointProperties,
157 typename vinvValueType, class ...vinvProperties>
158 void
159 Basis_HVOL_TRI_Cn_FEM::
160 getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
161  const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
162  const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
163  const EOperator operatorType) {
164  typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
165  typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
166  typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
167  typedef typename ExecSpace<typename inputPointViewType::execution_space,SpT>::ExecSpaceType ExecSpaceType;
168 
169  // loopSize corresponds to cardinality
170  const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
171  const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
172  const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
173  Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
174 
175  typedef typename inputPointViewType::value_type inputPointType;
176 
177  const ordinal_type cardinality = outputValues.extent(0);
178  const ordinal_type spaceDim = 2;
179 
180  auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
181  typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
182 
183  switch (operatorType) {
184  case OPERATOR_VALUE: {
185  workViewType work(Kokkos::view_alloc("Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
186  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
187  OPERATOR_VALUE,numPtsPerEval> FunctorType;
188  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
189  break;
190  }
191  case OPERATOR_GRAD:
192  case OPERATOR_D1: {
193  workViewType work(Kokkos::view_alloc("Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
194  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
195  OPERATOR_D1,numPtsPerEval> FunctorType;
196  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
197  break;
198  }
199  case OPERATOR_D2: {
200  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
201  OPERATOR_D2,numPtsPerEval> FunctorType;
202  workViewType work(Kokkos::view_alloc("Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
203  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
204  break;
205  }
206  /* case OPERATOR_D3: {
207  typedef Functor<outputValueViewType,inputPointViewType,vinvViewType,
208  OPERATOR_D3,numPtsPerEval> FunctorType;
209  workViewType work(Kokkos::view_alloc("Basis_HVOL_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0), outputValues.extent(2));
210  Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
211  break;
212  }*/
213  default: {
214  INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
215  ">>> ERROR (Basis_HVOL_TRI_Cn_FEM): Operator type not implemented" );
216  }
217  }
218 }
219 }
220 
221 // -------------------------------------------------------------------------------------
222 template<typename SpT, typename OT, typename PT>
224 Basis_HVOL_TRI_Cn_FEM( const ordinal_type order,
225  const EPointType pointType ) {
226 
227  constexpr ordinal_type spaceDim = 2;
228 
229  this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order); // bigN
230  this->basisDegree_ = order; // small n
231  this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
232  this->basisType_ = BASIS_FEM_FIAT;
233  this->basisCoordinates_ = COORDINATES_CARTESIAN;
234  this->functionSpace_ = FUNCTION_SPACE_HVOL;
235 
236  const ordinal_type card = this->basisCardinality_;
237 
238  // points are computed in the host and will be copied
239  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
240  dofCoords("HVOL::Tri::Cn::dofCoords", card, spaceDim);
241 
242  // construct lattice (only internal nodes for HVOL element)
243  const ordinal_type offset = 1;
244  PointTools::getLattice( dofCoords,
245  this->basisCellTopology_,
246  order+spaceDim+offset, offset,
247  pointType );
248 
249  this->dofCoords_ = Kokkos::create_mirror_view(typename SpT::memory_space(), dofCoords);
250  Kokkos::deep_copy(this->dofCoords_, dofCoords);
251 
252  // form Vandermonde matrix. Actually, this is the transpose of the VDM,
253  // so we transpose on copy below.
254  const ordinal_type lwork = card*card;
255  Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
256  vmat("HVOL::Tri::Cn::vmat", card, card),
257  work("HVOL::Tri::Cn::work", lwork),
258  ipiv("HVOL::Tri::Cn::ipiv", card);
259 
260  Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(vmat, dofCoords, order, OPERATOR_VALUE);
261 
262  ordinal_type info = 0;
263  Teuchos::LAPACK<ordinal_type,scalarType> lapack;
264 
265  lapack.GETRF(card, card,
266  vmat.data(), vmat.stride_1(),
267  (ordinal_type*)ipiv.data(),
268  &info);
269 
270  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
271  std::runtime_error ,
272  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
273 
274  lapack.GETRI(card,
275  vmat.data(), vmat.stride_1(),
276  (ordinal_type*)ipiv.data(),
277  work.data(), lwork,
278  &info);
279 
280  INTREPID2_TEST_FOR_EXCEPTION( info != 0,
281  std::runtime_error ,
282  ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
283 
284  // create host mirror
285  Kokkos::DynRankView<scalarType,typename SpT::array_layout,Kokkos::HostSpace>
286  vinv("HVOL::Line::Cn::vinv", card, card);
287 
288  for (ordinal_type i=0;i<card;++i)
289  for (ordinal_type j=0;j<card;++j)
290  vinv(i,j) = vmat(j,i);
291 
292  this->vinv_ = Kokkos::create_mirror_view(typename SpT::memory_space(), vinv);
293  Kokkos::deep_copy(this->vinv_ , vinv);
294 
295  // initialize tags
296  {
297  // Basis-dependent initializations
298  constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
299  const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
300  const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
301  const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
302 
303  constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
304  ordinal_type tags[maxCard][tagSize];
305 
306  const ordinal_type
307  numElemDof = this->basisCardinality_; //all the degrees of freedom are internal.
308 
309  ordinal_type elemId = 0;
310  for (ordinal_type i=0;i<this->basisCardinality_;++i) {
311  // elem
312  tags[i][0] = spaceDim; // intr dof
313  tags[i][1] = 0; // intr id
314  tags[i][2] = elemId++; // local dof id
315  tags[i][3] = numElemDof; // total vert dof
316  }
317 
318  OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
319 
320  // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
321  // tags are constructed on host
322  this->setOrdinalTagData(this->tagToOrdinal_,
323  this->ordinalToTag_,
324  tagView,
325  this->basisCardinality_,
326  tagSize,
327  posScDim,
328  posScOrd,
329  posDfOrd);
330  }
331 }
332 } // namespace Intrepid2
333 #endif
Kokkos::View< ordinal_type *, typename ExecSpaceType::array_layout, Kokkos::HostSpace > OrdinalTypeArray1DHost
View type for 1d host array.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
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 (currently disabled for other ce...
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.