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
53namespace Intrepid2 {
54
55// -------------------------------------------------------------------------------------
56namespace Impl {
57
58template<EOperator opType>
59template<typename OutputViewType,
60typename inputViewType,
61typename workViewType,
62typename vinvViewType>
63KOKKOS_INLINE_FUNCTION
64void
65Basis_HVOL_TRI_Cn_FEM::Serial<opType>::
66getValues( 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
154template<typename DT, ordinal_type numPtsPerEval,
155typename outputValueValueType, class ...outputValueProperties,
156typename inputPointValueType, class ...inputPointProperties,
157typename vinvValueType, class ...vinvProperties>
158void
159Basis_HVOL_TRI_Cn_FEM::
160getValues( 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,typename DT::execution_space>::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// -------------------------------------------------------------------------------------
222template<typename DT, typename OT, typename PT>
224Basis_HVOL_TRI_Cn_FEM( const ordinal_type order,
225 const EPointType pointType ) {
226
227 constexpr ordinal_type spaceDim = 2;
228
229 this->pointType_ = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
230 this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order); // bigN
231 this->basisDegree_ = order; // small n
232 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
233 this->basisType_ = BASIS_FEM_LAGRANGIAN;
234 this->basisCoordinates_ = COORDINATES_CARTESIAN;
235 this->functionSpace_ = FUNCTION_SPACE_HVOL;
236
237 const ordinal_type card = this->basisCardinality_;
238
239 // points are computed in the host and will be copied
240 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
241 dofCoords("HVOL::Tri::Cn::dofCoords", card, spaceDim);
242
243 // construct lattice (only internal nodes for HVOL element)
244 const ordinal_type offset = 1;
245 PointTools::getLattice( dofCoords,
246 this->basisCellTopology_,
247 order+spaceDim+offset, offset,
248 this->pointType_ );
249
250 this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
251 Kokkos::deep_copy(this->dofCoords_, dofCoords);
252
253 // form Vandermonde matrix. Actually, this is the transpose of the VDM,
254 // so we transpose on copy below.
255 const ordinal_type lwork = card*card;
256 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
257 vmat("HVOL::Tri::Cn::vmat", card, card),
258 work("HVOL::Tri::Cn::work", lwork),
259 ipiv("HVOL::Tri::Cn::ipiv", card);
260
261 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(vmat, dofCoords, order, OPERATOR_VALUE);
262
263 ordinal_type info = 0;
264 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
265
266 lapack.GETRF(card, card,
267 vmat.data(), vmat.stride_1(),
268 (ordinal_type*)ipiv.data(),
269 &info);
270
271 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
272 std::runtime_error ,
273 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
274
275 lapack.GETRI(card,
276 vmat.data(), vmat.stride_1(),
277 (ordinal_type*)ipiv.data(),
278 work.data(), lwork,
279 &info);
280
281 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
282 std::runtime_error ,
283 ">>> ERROR: (Intrepid2::Basis_HVOL_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
284
285 // create host mirror
286 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
287 vinv("HVOL::Line::Cn::vinv", card, card);
288
289 for (ordinal_type i=0;i<card;++i)
290 for (ordinal_type j=0;j<card;++j)
291 vinv(i,j) = vmat(j,i);
292
293 this->vinv_ = Kokkos::create_mirror_view(typename DT::memory_space(), vinv);
294 Kokkos::deep_copy(this->vinv_ , vinv);
295
296 // initialize tags
297 {
298 // Basis-dependent initializations
299 constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
300 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
301 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
302 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
303
304 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
305 ordinal_type tags[maxCard][tagSize];
306
307 const ordinal_type
308 numElemDof = this->basisCardinality_; //all the degrees of freedom are internal.
309
310 ordinal_type elemId = 0;
311 for (ordinal_type i=0;i<this->basisCardinality_;++i) {
312 // elem
313 tags[i][0] = spaceDim; // intr dof
314 tags[i][1] = 0; // intr id
315 tags[i][2] = elemId++; // local dof id
316 tags[i][3] = numElemDof; // total vert dof
317 }
318
319 OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
320
321 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
322 // tags are constructed on host
323 this->setOrdinalTagData(this->tagToOrdinal_,
324 this->ordinalToTag_,
325 tagView,
326 this->basisCardinality_,
327 tagSize,
328 posScDim,
329 posScOrd,
330 posDfOrd);
331 }
332}
333} // namespace Intrepid2
334#endif
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Basis_HVOL_TRI_Cn_FEM(const ordinal_type order, const EPointType pointType=POINTTYPE_EQUISPACED)
Constructor.
static constexpr ordinal_type MaxOrder
The maximum reconstruction order.
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...