Intrepid2
Intrepid2_HGRAD_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
49#ifndef __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
50#define __INTREPID2_HGRAD_TRI_CN_FEM_DEF_HPP__
51
53
54namespace Intrepid2 {
55
56// -------------------------------------------------------------------------------------
57namespace Impl {
58
59template<EOperator opType>
60template<typename OutputViewType,
61typename inputViewType,
62typename workViewType,
63typename vinvViewType>
64KOKKOS_INLINE_FUNCTION
65void
66Basis_HGRAD_TRI_Cn_FEM::Serial<opType>::
67getValues( OutputViewType output,
68 const inputViewType input,
69 workViewType work,
70 const vinvViewType vinv ) {
71
72 constexpr ordinal_type spaceDim = 2;
73 const ordinal_type
74 card = vinv.extent(0),
75 npts = input.extent(0);
76
77 // compute order
78 ordinal_type order = 0;
79 for (ordinal_type p=0;p<=Parameters::MaxOrder;++p) {
80 if (card == Intrepid2::getPnCardinality<spaceDim>(p) ) {
81 order = p;
82 break;
83 }
84 }
85
86 typedef typename Kokkos::DynRankView<typename workViewType::value_type, typename workViewType::memory_space> viewType;
87 auto vcprop = Kokkos::common_view_alloc_prop(work);
88 auto ptr = work.data();
89
90 switch (opType) {
91 case OPERATOR_VALUE: {
92 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts);
93 viewType dummyView;
94
95 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
96 Serial<opType>::getValues(phis, input, dummyView, order);
97
98 for (ordinal_type i=0;i<card;++i)
99 for (ordinal_type j=0;j<npts;++j) {
100 output.access(i,j) = 0.0;
101 for (ordinal_type k=0;k<card;++k)
102 output.access(i,j) += vinv(k,i)*phis.access(k,j);
103 }
104 break;
105 }
106 case OPERATOR_GRAD:
107 case OPERATOR_D1: {
108 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
109 ptr += card*npts*spaceDim*get_dimension_scalar(work);
110 const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
111 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
112 Serial<opType>::getValues(phis, input, workView, order);
113
114 for (ordinal_type i=0;i<card;++i)
115 for (ordinal_type j=0;j<npts;++j)
116 for (ordinal_type k=0;k<spaceDim;++k) {
117 output.access(i,j,k) = 0.0;
118 for (ordinal_type l=0;l<card;++l)
119 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
120 }
121 break;
122 }
123 case OPERATOR_D2:
124 case OPERATOR_D3:
125 case OPERATOR_D4:
126 case OPERATOR_D5:
127 case OPERATOR_D6:
128 case OPERATOR_D7:
129 case OPERATOR_D8:
130 case OPERATOR_D9:
131 case OPERATOR_D10: {
132 const ordinal_type dkcard = getDkCardinality<opType,spaceDim>(); //(orDn + 1);
133 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, dkcard);
134 viewType dummyView;
135
136 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
137 Serial<opType>::getValues(phis, input, dummyView, order);
138
139 for (ordinal_type i=0;i<card;++i)
140 for (ordinal_type j=0;j<npts;++j)
141 for (ordinal_type k=0;k<dkcard;++k) {
142 output.access(i,j,k) = 0.0;
143 for (ordinal_type l=0;l<card;++l)
144 output.access(i,j,k) += vinv(l,i)*phis.access(l,j,k);
145 }
146 break;
147 }
148 case OPERATOR_CURL: { // only works in 2d. first component is -d/dy, second is d/dx
149 const viewType phis(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim);
150 ptr += card*npts*spaceDim*get_dimension_scalar(work);
151 const viewType workView(Kokkos::view_wrap(ptr, vcprop), card, npts, spaceDim+1);
152
153
154 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::
155 Serial<OPERATOR_D1>::getValues(phis, input, workView, order);
156
157 for (ordinal_type i=0;i<card;++i)
158 for (ordinal_type j=0;j<npts;++j) {
159 output.access(i,j,0) = 0.0;
160 for (ordinal_type l=0;l<card;++l)
161 output.access(i,j,0) += vinv(l,i)*phis.access(l,j,1);
162 output.access(i,j,1) = 0.0;
163 for (ordinal_type l=0;l<card;++l)
164 output.access(i,j,1) -= vinv(l,i)*phis.access(l,j,0);
165 }
166 break;
167 }
168 default: {
169 INTREPID2_TEST_FOR_ABORT( true,
170 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented");
171 }
172 }
173}
174
175template<typename DT, ordinal_type numPtsPerEval,
176typename outputValueValueType, class ...outputValueProperties,
177typename inputPointValueType, class ...inputPointProperties,
178typename vinvValueType, class ...vinvProperties>
179void
180Basis_HGRAD_TRI_Cn_FEM::
181getValues( Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValues,
182 const Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPoints,
183 const Kokkos::DynRankView<vinvValueType, vinvProperties...> vinv,
184 const EOperator operatorType) {
185 typedef Kokkos::DynRankView<outputValueValueType,outputValueProperties...> outputValueViewType;
186 typedef Kokkos::DynRankView<inputPointValueType, inputPointProperties...> inputPointViewType;
187 typedef Kokkos::DynRankView<vinvValueType, vinvProperties...> vinvViewType;
188 typedef typename ExecSpace<typename inputPointViewType::execution_space,typename DT::execution_space>::ExecSpaceType ExecSpaceType;
189
190 // loopSize corresponds to cardinality
191 const auto loopSizeTmp1 = (inputPoints.extent(0)/numPtsPerEval);
192 const auto loopSizeTmp2 = (inputPoints.extent(0)%numPtsPerEval != 0);
193 const auto loopSize = loopSizeTmp1 + loopSizeTmp2;
194 Kokkos::RangePolicy<ExecSpaceType,Kokkos::Schedule<Kokkos::Static> > policy(0, loopSize);
195
196 typedef typename inputPointViewType::value_type inputPointType;
197
198 const ordinal_type cardinality = outputValues.extent(0);
199 const ordinal_type spaceDim = 2;
200
201 auto vcprop = Kokkos::common_view_alloc_prop(inputPoints);
202 typedef typename Kokkos::DynRankView< inputPointType, typename inputPointViewType::memory_space> workViewType;
203
204 switch (operatorType) {
205 case OPERATOR_VALUE: {
206 workViewType work(Kokkos::view_alloc("Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0));
207 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
208 OPERATOR_VALUE,numPtsPerEval> FunctorType;
209 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
210 break;
211 }
212 case OPERATOR_GRAD:
213 case OPERATOR_D1: {
214 workViewType work(Kokkos::view_alloc("Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
215 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
216 OPERATOR_D1,numPtsPerEval> FunctorType;
217 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
218 break;
219 }
220 case OPERATOR_CURL: {
221 workViewType work(Kokkos::view_alloc("Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*(2*spaceDim+1), inputPoints.extent(0));
222 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
223 OPERATOR_CURL,numPtsPerEval> FunctorType;
224 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
225 break;
226 }
227 case OPERATOR_D2: {
228 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType,
229 OPERATOR_D2,numPtsPerEval> FunctorType;
230 workViewType work(Kokkos::view_alloc("Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality*outputValues.extent(2), inputPoints.extent(0));
231 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
232 break;
233 }
234 /* case OPERATOR_D3: {
235 typedef Functor<outputValueViewType,inputPointViewType,vinvViewType, workViewType
236 OPERATOR_D3,numPtsPerEval> FunctorType;
237 workViewType work(Kokkos::view_alloc("Basis_HGRAD_TRI_Cn_FEM::getValues::work", vcprop), cardinality, inputPoints.extent(0), outputValues.extent(2));
238 Kokkos::parallel_for( policy, FunctorType(outputValues, inputPoints, vinv, work) );
239 break;
240 }*/
241 default: {
242 INTREPID2_TEST_FOR_EXCEPTION( true , std::invalid_argument,
243 ">>> ERROR (Basis_HGRAD_TRI_Cn_FEM): Operator type not implemented" );
244 }
245 }
246}
247}
248
249// -------------------------------------------------------------------------------------
250template<typename DT, typename OT, typename PT>
252Basis_HGRAD_TRI_Cn_FEM( const ordinal_type order,
253 const EPointType pointType ) {
254 constexpr ordinal_type spaceDim = 2;
255
256 this->basisCardinality_ = Intrepid2::getPnCardinality<spaceDim>(order); // bigN
257 this->basisDegree_ = order; // small n
258 this->basisCellTopology_ = shards::CellTopology(shards::getCellTopologyData<shards::Triangle<3> >() );
259 this->basisType_ = BASIS_FEM_LAGRANGIAN;
260 this->basisCoordinates_ = COORDINATES_CARTESIAN;
261 this->functionSpace_ = FUNCTION_SPACE_HGRAD;
262
263 pointType_ = (pointType == POINTTYPE_DEFAULT) ? POINTTYPE_EQUISPACED : pointType;
264 const ordinal_type card = this->basisCardinality_;
265
266 // points are computed in the host and will be copied
267 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
268 dofCoords("Hgrad::Tri::Cn::dofCoords", card, spaceDim);
269
270 // construct lattice
271 const ordinal_type offset = 0;
272 PointTools::getLattice( dofCoords,
273 this->basisCellTopology_,
274 order, offset,
275 pointType_ );
276
277 this->dofCoords_ = Kokkos::create_mirror_view(typename DT::memory_space(), dofCoords);
278 Kokkos::deep_copy(this->dofCoords_, dofCoords);
279
280 // form Vandermonde matrix. Actually, this is the transpose of the VDM,
281 // so we transpose on copy below.
282 const ordinal_type lwork = card*card;
283 Kokkos::DynRankView<scalarType,Kokkos::LayoutLeft,Kokkos::HostSpace>
284 vmat("Hgrad::Tri::Cn::vmat", card, card),
285 work("Hgrad::Tri::Cn::work", lwork),
286 ipiv("Hgrad::Tri::Cn::ipiv", card);
287
288 Impl::Basis_HGRAD_TRI_Cn_FEM_ORTH::getValues<Kokkos::HostSpace::execution_space,Parameters::MaxNumPtsPerBasisEval>(vmat, dofCoords, order, OPERATOR_VALUE);
289
290 ordinal_type info = 0;
291 Teuchos::LAPACK<ordinal_type,scalarType> lapack;
292
293 lapack.GETRF(card, card,
294 vmat.data(), vmat.stride_1(),
295 (ordinal_type*)ipiv.data(),
296 &info);
297
298 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
299 std::runtime_error ,
300 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRF returns nonzero info." );
301
302 lapack.GETRI(card,
303 vmat.data(), vmat.stride_1(),
304 (ordinal_type*)ipiv.data(),
305 work.data(), lwork,
306 &info);
307
308 INTREPID2_TEST_FOR_EXCEPTION( info != 0,
309 std::runtime_error ,
310 ">>> ERROR: (Intrepid2::Basis_HGRAD_TRI_Cn_FEM) lapack.GETRI returns nonzero info." );
311
312 // create host mirror
313 Kokkos::DynRankView<scalarType,typename DT::execution_space::array_layout,Kokkos::HostSpace>
314 vinv("Hgrad::Line::Cn::vinv", card, card);
315
316 for (ordinal_type i=0;i<card;++i)
317 for (ordinal_type j=0;j<card;++j)
318 vinv(i,j) = vmat(j,i);
319
320 this->vinv_ = Kokkos::create_mirror_view(typename DT::memory_space(), vinv);
321 Kokkos::deep_copy(this->vinv_ , vinv);
322
323 // initialize tags
324 {
325 // Basis-dependent initializations
326 constexpr ordinal_type tagSize = 4; // size of DoF tag, i.e., number of fields in the tag
327 const ordinal_type posScDim = 0; // position in the tag, counting from 0, of the subcell dim
328 const ordinal_type posScOrd = 1; // position in the tag, counting from 0, of the subcell ordinal
329 const ordinal_type posDfOrd = 2; // position in the tag, counting from 0, of DoF ordinal relative to the subcell
330
331 // 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.
332 // 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.)
333 INTREPID2_TEST_FOR_EXCEPTION( order > Parameters::MaxOrder, std::invalid_argument, "polynomial order exceeds the max supported by this class");
334 constexpr ordinal_type maxCard = Intrepid2::getPnCardinality<spaceDim, Parameters::MaxOrder>();
335 ordinal_type tags[maxCard][tagSize];
336
337 const ordinal_type
338 numEdgeDof = Intrepid2::getPnCardinality<1>(order-2),
339 numElemDof = (order > 2 ? Intrepid2::getPnCardinality<2>(order-3) : 0);
340
341 scalarType xi0, xi1, xi2;
342 const scalarType eps = threshold();
343
344 ordinal_type edgeId[3] = {}, elemId = 0;
345 for (ordinal_type i=0;i<card;++i) {
346
347 // compute barycentric coordinates
348 const auto x = dofCoords(i,0);
349 const auto y = dofCoords(i,1);
350 xi0 = 1.0 - x - y;
351 xi1= x;
352 xi2= y;
353
354 // vertex
355 if ((1.0 - xi0) < eps) { // vert 0
356 tags[i][0] = 0; // vertex dof
357 tags[i][1] = 0; // vertex id
358 tags[i][2] = 0; // local dof id
359 tags[i][3] = 1; // total vert dof
360 }
361 else if ((1.0 - xi1) < eps) { // vert 1
362 tags[i][0] = 0; // vertex dof
363 tags[i][1] = 1; // vertex id
364 tags[i][2] = 0; // local dof id
365 tags[i][3] = 1; // total vert dof
366 }
367 else if ((1.0 - xi2) < eps) { // vert 2
368 tags[i][0] = 0; // vertex dof
369 tags[i][1] = 2; // vertex id
370 tags[i][2] = 0; // local dof id
371 tags[i][3] = 1; // total vert dof
372 }
373 else if (xi2 < eps) { // edge 0
374 tags[i][0] = 1; // edge dof
375 tags[i][1] = 0; // edge id
376 tags[i][2] = edgeId[0]++; // local dof id
377 tags[i][3] = numEdgeDof; // total vert dof
378 }
379 else if (xi0 < eps) { // edge 1
380 tags[i][0] = 1; // edge dof
381 tags[i][1] = 1; // edge id
382 tags[i][2] = edgeId[1]++; // local dof id
383 tags[i][3] = numEdgeDof; // total vert dof
384 }
385 else if (xi1 < eps) { // edge 2
386 tags[i][0] = 1; // edge dof
387 tags[i][1] = 2; // edge id
388 tags[i][2] = edgeId[2]++; // local dof id
389 tags[i][3] = numEdgeDof; // total vert dof
390 }
391 else { // elem
392 tags[i][0] = 2; // intr dof
393 tags[i][1] = 0; // intr id
394 tags[i][2] = elemId++; // local dof id
395 tags[i][3] = numElemDof; // total vert dof
396 }
397 }
398
399 OrdinalTypeArray1DHost tagView(&tags[0][0], card*tagSize);
400
401 // Basis-independent function sets tag and enum data in tagToOrdinal_ and ordinalToTag_ arrays:
402 // tags are constructed on host
403 this->setOrdinalTagData(this->tagToOrdinal_,
404 this->ordinalToTag_,
405 tagView,
406 this->basisCardinality_,
407 tagSize,
408 posScDim,
409 posScOrd,
410 posDfOrd);
411 }
412}
413} // namespace Intrepid2
414#endif
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM_ORTH class.
Basis_HGRAD_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...