Intrepid2
Intrepid2_ProjectionToolsDefHCURL.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_PROJECTIONTOOLSDEFHCURL_HPP__
50 #define __INTREPID2_PROJECTIONTOOLSDEFHCURL_HPP__
51 
53 #include "Intrepid2_ArrayTools.hpp"
55 
56 
57 namespace Intrepid2 {
58 namespace Experimental {
59 
60 template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
62  const ViewType1 basisTanAtBasisEPoints_;
63  const ViewType1 basisAtBasisEPoints_;
64  const ViewType2 basisEWeights_;
65  const ViewType1 wTanBasisAtBasisEPoints_;
66  const ViewType2 targetEWeights_;
67  const ViewType1 basisAtTargetEPoints_;
68  const ViewType1 wTanBasisAtTargetEPoints_;
69  const ViewType3 tagToOrdinal_;
70  const ViewType4 targetAtTargetEPoints_;
71  const ViewType1 targetTanAtTargetEPoints_;
72  const ViewType1 refEdgesTangent_;
73  ordinal_type edgeCardinality_;
74  ordinal_type offsetBasis_;
75  ordinal_type offsetTarget_;
76  ordinal_type edgeDim_;
77  ordinal_type dim_;
78  ordinal_type iedge_;
79 
80  ComputeBasisCoeffsOnEdges_HCurl(const ViewType1 basisTanAtBasisEPoints,
81  const ViewType1 basisAtBasisEPoints, const ViewType2 basisEWeights, const ViewType1 wTanBasisAtBasisEPoints, const ViewType2 targetEWeights,
82  const ViewType1 basisAtTargetEPoints, const ViewType1 wTanBasisAtTargetEPoints, const ViewType3 tagToOrdinal,
83  const ViewType4 targetAtTargetEPoints, const ViewType1 targetTanAtTargetEPoints,
84  const ViewType1 refEdgesTangent, ordinal_type edgeCardinality, ordinal_type offsetBasis,
85  ordinal_type offsetTarget, ordinal_type edgeDim,
86  ordinal_type dim, ordinal_type iedge) :
87  basisTanAtBasisEPoints_(basisTanAtBasisEPoints),
88  basisAtBasisEPoints_(basisAtBasisEPoints), basisEWeights_(basisEWeights), wTanBasisAtBasisEPoints_(wTanBasisAtBasisEPoints), targetEWeights_(targetEWeights),
89  basisAtTargetEPoints_(basisAtTargetEPoints), wTanBasisAtTargetEPoints_(wTanBasisAtTargetEPoints),
90  tagToOrdinal_(tagToOrdinal), targetAtTargetEPoints_(targetAtTargetEPoints),
91  targetTanAtTargetEPoints_(targetTanAtTargetEPoints),
92  refEdgesTangent_(refEdgesTangent), edgeCardinality_(edgeCardinality), offsetBasis_(offsetBasis),
93  offsetTarget_(offsetTarget), edgeDim_(edgeDim), dim_(dim), iedge_(iedge)
94  {}
95 
96  void
97  KOKKOS_INLINE_FUNCTION
98  operator()(const ordinal_type ic) const {
99 
100  ordinal_type numBasisEPoints = basisEWeights_.extent(0);
101  ordinal_type numTargetEPoints = targetEWeights_.extent(0);
102  for(ordinal_type j=0; j <edgeCardinality_; ++j) {
103  ordinal_type jdof = tagToOrdinal_(edgeDim_, iedge_, j);
104  for(ordinal_type iq=0; iq <numBasisEPoints; ++iq) {
105  for(ordinal_type d=0; d <dim_; ++d)
106  basisTanAtBasisEPoints_(ic,j,iq) += refEdgesTangent_(iedge_,d)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
107  wTanBasisAtBasisEPoints_(ic,j,iq) = basisTanAtBasisEPoints_(ic,j,iq)*basisEWeights_(iq);
108  }
109 
110  for(ordinal_type iq=0; iq <numTargetEPoints; ++iq) {
111  typename ViewType2::value_type tmp = 0;
112  for(ordinal_type d=0; d <dim_; ++d)
113  tmp += refEdgesTangent_(iedge_,d)*basisAtTargetEPoints_(ic,jdof,offsetTarget_+iq,d);
114  wTanBasisAtTargetEPoints_(ic,j,iq) = tmp*targetEWeights_(iq);
115  }
116  }
117  for(ordinal_type iq=0; iq <numTargetEPoints; ++iq)
118  for(ordinal_type d=0; d <dim_; ++d)
119  targetTanAtTargetEPoints_(ic,iq) += refEdgesTangent_(iedge_,d)*targetAtTargetEPoints_(ic,offsetTarget_+iq,d);
120  }
121 };
122 
123 
124 template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4,
125 typename ViewType5, typename ViewType6, typename ViewType7, typename ViewType8>
127  const ViewType1 basisCoeffs_;
128  const ViewType2 orts_;
129  const ViewType3 negPartialProjTan_;
130  const ViewType3 negPartialProjCurlNormal_;
131  const ViewType3 hgradBasisGradAtBasisEPoints_;
132  const ViewType3 wHgradBasisGradAtBasisEPoints_;
133  const ViewType3 basisCurlAtBasisCurlEPoints_;
134  const ViewType3 basisCurlNormalAtBasisCurlEPoints_;
135  const ViewType3 basisAtBasisEPoints_;
136  const ViewType3 normalTargetCurlAtTargetEPoints_;
137  const ViewType3 basisTanAtBasisEPoints_;
138  const ViewType3 hgradBasisGradAtTargetEPoints_;
139  const ViewType3 wHgradBasisGradAtTargetEPoints_;
140  const ViewType3 wNormalBasisCurlAtBasisCurlEPoints_;
141  const ViewType3 basisCurlAtTargetCurlEPoints_;
142  const ViewType3 wNormalBasisCurlBasisAtTargetCurlEPoints_;
143  const ViewType4 targetAtTargetEPoints_;
144  const ViewType3 targetTanAtTargetEPoints_;
145  const ViewType4 targetCurlAtTargetCurlEPoints_;
146  const ViewType5 basisEWeights_;
147  const ViewType5 targetEWeights_;
148  const ViewType5 basisCurlEWeights_;
149  const ViewType5 targetCurlEWeights_;
150  const ViewType6 tagToOrdinal_;
151  const ViewType6 hGradTagToOrdinal_;
152  const ViewType7 refTopologyKey_;
153  const ViewType3 refFacesNormal_;
154  const ViewType3 refFacesTangents_;
155  const ViewType8 computedDofs_;
156  ordinal_type offsetBasis_;
157  ordinal_type offsetBasisCurl_;
158  ordinal_type offsetTarget_;
159  ordinal_type offsetTargetCurl_;
160  ordinal_type iface_;
161  ordinal_type hgradCardinality_;
162  ordinal_type numFaces_;
163  ordinal_type numFaceDofs_;
164  ordinal_type numEdgeDofs_;
165  ordinal_type faceDim_;
166  ordinal_type dim_;
167 
168 
169 
170 
171  ComputeBasisCoeffsOnFaces_HCurl(const ViewType1 basisCoeffs,
172  const ViewType2 orts, const ViewType3 negPartialProjTan, const ViewType3 negPartialProjCurlNormal,
173  const ViewType3 hgradBasisGradAtBasisEPoints, const ViewType3 wHgradBasisGradAtBasisEPoints,
174  const ViewType3 basisCurlAtBasisCurlEPoints, const ViewType3 basisCurlNormalAtBasisCurlEPoints,
175  const ViewType3 basisAtBasisEPoints,
176  const ViewType3 normalTargetCurlAtTargetEPoints,
177  const ViewType3 basisTanAtBasisEPoints,
178  const ViewType3 hgradBasisGradAtTargetEPoints, const ViewType3 wHgradBasisGradAtTargetEPoints,
179  const ViewType3 wNormalBasisCurlAtBasisCurlEPoints, const ViewType3 basisCurlAtTargetCurlEPoints,
180  const ViewType3 wNormalBasisCurlBasisAtTargetCurlEPoints, const ViewType4 targetAtTargetEPoints,
181  const ViewType3 targetTanAtTargetEPoints, const ViewType4 targetCurlAtTargetCurlEPoints,
182  const ViewType5 basisEWeights, const ViewType5 targetEWeights,
183  const ViewType5 basisCurlEWeights, const ViewType5 targetCurlEWeights, const ViewType6 tagToOrdinal,
184  const ViewType6 hGradTagToOrdinal, const ViewType7 refTopologyKey,
185  const ViewType3 refFacesNormal, const ViewType3 refFacesTangents,
186  const ViewType8 computedDofs, ordinal_type offsetBasis,
187  ordinal_type offsetBasisCurl, ordinal_type offsetTarget,
188  ordinal_type offsetTargetCurl, ordinal_type iface,
189  ordinal_type hgradCardinality, ordinal_type numFaces,
190  ordinal_type numFaceDofs, ordinal_type numEdgeDofs,
191  ordinal_type faceDim, ordinal_type dim):
192  basisCoeffs_(basisCoeffs),
193  orts_(orts), negPartialProjTan_(negPartialProjTan), negPartialProjCurlNormal_(negPartialProjCurlNormal),
194  hgradBasisGradAtBasisEPoints_(hgradBasisGradAtBasisEPoints), wHgradBasisGradAtBasisEPoints_(wHgradBasisGradAtBasisEPoints),
195  basisCurlAtBasisCurlEPoints_(basisCurlAtBasisCurlEPoints), basisCurlNormalAtBasisCurlEPoints_(basisCurlNormalAtBasisCurlEPoints),
196  basisAtBasisEPoints_(basisAtBasisEPoints),
197  normalTargetCurlAtTargetEPoints_(normalTargetCurlAtTargetEPoints), basisTanAtBasisEPoints_(basisTanAtBasisEPoints),
198  hgradBasisGradAtTargetEPoints_(hgradBasisGradAtTargetEPoints), wHgradBasisGradAtTargetEPoints_(wHgradBasisGradAtTargetEPoints),
199  wNormalBasisCurlAtBasisCurlEPoints_(wNormalBasisCurlAtBasisCurlEPoints), basisCurlAtTargetCurlEPoints_(basisCurlAtTargetCurlEPoints),
200  wNormalBasisCurlBasisAtTargetCurlEPoints_(wNormalBasisCurlBasisAtTargetCurlEPoints), targetAtTargetEPoints_(targetAtTargetEPoints),
201  targetTanAtTargetEPoints_(targetTanAtTargetEPoints), targetCurlAtTargetCurlEPoints_(targetCurlAtTargetCurlEPoints),
202  basisEWeights_(basisEWeights), targetEWeights_(targetEWeights),
203  basisCurlEWeights_(basisCurlEWeights), targetCurlEWeights_(targetCurlEWeights), tagToOrdinal_(tagToOrdinal),
204  hGradTagToOrdinal_(hGradTagToOrdinal), refTopologyKey_(refTopologyKey),
205  refFacesNormal_(refFacesNormal), refFacesTangents_(refFacesTangents),
206  computedDofs_(computedDofs), offsetBasis_(offsetBasis),
207  offsetBasisCurl_(offsetBasisCurl), offsetTarget_(offsetTarget),
208  offsetTargetCurl_(offsetTargetCurl), iface_(iface),
209  hgradCardinality_(hgradCardinality), numFaces_(numFaces),
210  numFaceDofs_(numFaceDofs), numEdgeDofs_(numEdgeDofs),
211  faceDim_(faceDim), dim_(dim){}
212 
213  void
214  KOKKOS_INLINE_FUNCTION
215  operator()(const ordinal_type ic) const {
216 
217  ordinal_type fOrt[6];
218  orts_(ic).getFaceOrientation(fOrt, numFaces_);
219 
220  ordinal_type ort = fOrt[iface_];
221  typename ViewType3::value_type ortJacData[4]; //faceDim x faceDim
222  auto ortJac = ViewType3(ortJacData, faceDim_, faceDim_);
223  Impl::OrientationTools::getJacobianOfOrientationMap(ortJac, refTopologyKey_(faceDim_,iface_), ort);
224 
225  ordinal_type numBasisEPoints = basisEWeights_.extent(0);
226  ordinal_type numTargetEPoints = targetEWeights_.extent(0);
227  for(ordinal_type j=0; j <hgradCardinality_; ++j) {
228  ordinal_type face_dof = hGradTagToOrdinal_(faceDim_, 0, j);
229  for(ordinal_type d=0; d <faceDim_; ++d) {
230  for(ordinal_type iq=0; iq <numBasisEPoints; ++iq)
231  wHgradBasisGradAtBasisEPoints_(ic, j, iq, d) = hgradBasisGradAtBasisEPoints_(face_dof,iq,d) * basisEWeights_(iq);
232  for(ordinal_type iq=0; iq <numTargetEPoints; ++iq)
233  wHgradBasisGradAtTargetEPoints_(ic,j,iq,d)= hgradBasisGradAtTargetEPoints_(face_dof,iq,d) * targetEWeights_(iq);
234  }
235  }
236 
237  //Note: we are not considering the jacobian of the orientation map for normals since it is simply a scalar term for the integrals and it does not affect the projection
238  ordinal_type numBasisCurlEPoints = basisCurlEWeights_.extent(0);
239  for(ordinal_type j=0; j <numFaceDofs_; ++j) {
240  ordinal_type jdof = tagToOrdinal_(faceDim_, iface_, j);
241  for(ordinal_type iq=0; iq <numBasisEPoints; ++iq) {
242  for(ordinal_type d=0; d <dim_; ++d) {
243  for(ordinal_type itan=0; itan <dim_-1; ++itan)
244  for(ordinal_type jtan=0; jtan <dim_-1; ++jtan)
245  basisTanAtBasisEPoints_(ic,j,iq,itan) += refFacesTangents_(iface_,d,jtan)*ortJac(jtan,itan)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
246  }
247  }
248  for(ordinal_type iq=0; iq <numBasisCurlEPoints; ++iq) {
249  for(ordinal_type d=0; d <dim_; ++d)
250  basisCurlNormalAtBasisCurlEPoints_(ic,j,iq) += refFacesNormal_(iface_,d)*basisCurlAtBasisCurlEPoints_(ic,jdof,offsetBasisCurl_+iq,d);
251  wNormalBasisCurlAtBasisCurlEPoints_(ic,j,iq) = basisCurlNormalAtBasisCurlEPoints_(ic,j,iq) * basisCurlEWeights_(iq);
252  }
253 
254  ordinal_type numTargetCurlEPoints = targetCurlEWeights_.extent(0);
255  for(ordinal_type iq=0; iq <numTargetCurlEPoints; ++iq) {
256  typename ViewType3::value_type tmp=0;
257  for(ordinal_type d=0; d <dim_; ++d)
258  tmp += refFacesNormal_(iface_,d)*basisCurlAtTargetCurlEPoints_(ic,jdof,offsetTargetCurl_+iq,d);
259  wNormalBasisCurlBasisAtTargetCurlEPoints_(ic,j,iq) = tmp*targetCurlEWeights_(iq);
260  }
261  }
262 
263  for(ordinal_type j=0; j <numEdgeDofs_; ++j) {
264  ordinal_type jdof = computedDofs_(j);
265  for(ordinal_type iq=0; iq <numBasisEPoints; ++iq)
266  for(ordinal_type d=0; d <dim_; ++d) {
267  negPartialProjCurlNormal_(ic,iq) -= refFacesNormal_(iface_,d)*basisCoeffs_(ic,jdof)*basisCurlAtBasisCurlEPoints_(ic,jdof,offsetBasisCurl_+iq,d);
268  for(ordinal_type itan=0; itan <dim_-1; ++itan)
269  for(ordinal_type jtan=0; jtan <dim_-1; ++jtan)
270  negPartialProjTan_(ic,iq,itan) -= refFacesTangents_(iface_,d,jtan)*ortJac(jtan,itan)*basisCoeffs_(ic,jdof)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
271  }
272  }
273 
274  ordinal_type numTargetCurlEPoints = targetCurlEWeights_.extent(0);
275  for(ordinal_type iq=0; iq <numTargetEPoints; ++iq)
276  for(ordinal_type d=0; d <dim_; ++d)
277  for(ordinal_type itan=0; itan <dim_-1; ++itan)
278  for(ordinal_type jtan=0; jtan <dim_-1; ++jtan)
279  targetTanAtTargetEPoints_(ic,iq,itan) += refFacesTangents_(iface_,d,jtan)*ortJac(jtan,itan)*targetAtTargetEPoints_(ic,offsetTarget_+iq,d);
280 
281  for(ordinal_type iq=0; iq <numTargetCurlEPoints; ++iq)
282  for(ordinal_type d=0; d <dim_; ++d)
283  normalTargetCurlAtTargetEPoints_(ic,iq) += refFacesNormal_(iface_,d)*targetCurlAtTargetCurlEPoints_(ic,offsetTargetCurl_+iq,d);
284  }
285 };
286 
287 
288 template<typename ViewType1, typename ViewType2, typename ViewType3,
289 typename ViewType4, typename ViewType5>
291  const ViewType1 basisCoeffs_;
292  const ViewType2 negPartialProj_;
293  const ViewType2 negPartialProjCurl_;
294  const ViewType2 cellBasisAtBasisEPoints_;
295  const ViewType2 cellBasisCurlAtBasisCurlEPoints_;
296  const ViewType2 basisAtBasisEPoints_;
297  const ViewType2 hgradBasisGradAtBasisEPoints_;
298  const ViewType2 basisCurlAtBasisCurlEPoints_;
299  const ViewType2 hgradBasisGradAtTargetEPoints_;
300  const ViewType2 basisCurlAtTargetCurlEPoints_;
301  const ViewType3 basisEWeights_;
302  const ViewType3 basisCurlEWeights_;
303  const ViewType2 wHgradBasisGradAtBasisEPoints_;
304  const ViewType2 wBasisCurlAtBasisCurlEPoints_;
305  const ViewType3 targetEWeights_;
306  const ViewType3 targetCurlEWeights_;
307  const ViewType2 wHgradBasisGradAtTargetEPoints_;
308  const ViewType2 wBasisCurlAtTargetCurlEPoints_;
309  const ViewType4 computedDofs_;
310  const ViewType5 tagToOrdinal_;
311  const ViewType5 hGradTagToOrdinal_;
312  ordinal_type numCellDofs_;
313  ordinal_type hgradCardinality_;
314  ordinal_type offsetBasis_;
315  ordinal_type offsetBasisCurl_;
316  ordinal_type offsetTargetCurl_;
317  ordinal_type numEdgeFaceDofs_;
318  ordinal_type dim_;
319  ordinal_type derDim_;
320 
321  ComputeBasisCoeffsOnCell_HCurl(const ViewType1 basisCoeffs, ViewType2 negPartialProj, ViewType2 negPartialProjCurl,
322  const ViewType2 cellBasisAtBasisEPoints, const ViewType2 cellBasisCurlAtBasisCurlEPoints,
323  const ViewType2 basisAtBasisEPoints, const ViewType2 hgradBasisGradAtBasisEPoints, const ViewType2 basisCurlAtBasisCurlEPoints,
324  const ViewType2 hgradBasisGradAtTargetEPoints, const ViewType2 basisCurlAtTargetCurlEPoints,
325  const ViewType3 basisEWeights, const ViewType3 basisCurlEWeights,
326  const ViewType2 wHgradBasisGradAtBasisEPoints, const ViewType2 wBasisCurlAtBasisCurlEPoints,
327  const ViewType3 targetEWeights, const ViewType3 targetCurlEWeights,
328  const ViewType2 wHgradBasisGradAtTargetEPoints,
329  const ViewType2 wBasisCurlAtTargetCurlEPoints, const ViewType4 computedDofs,
330  const ViewType5 tagToOrdinal, const ViewType5 hGradTagToOrdinal,
331  ordinal_type numCellDofs, ordinal_type hgradCardinality,
332  ordinal_type offsetBasis, ordinal_type offsetBasisCurl, ordinal_type offsetTargetCurl,
333  ordinal_type numEdgeFaceDofs, ordinal_type dim, ordinal_type derDim) :
334  basisCoeffs_(basisCoeffs), negPartialProj_(negPartialProj), negPartialProjCurl_(negPartialProjCurl),
335  cellBasisAtBasisEPoints_(cellBasisAtBasisEPoints), cellBasisCurlAtBasisCurlEPoints_(cellBasisCurlAtBasisCurlEPoints),
336  basisAtBasisEPoints_(basisAtBasisEPoints), hgradBasisGradAtBasisEPoints_(hgradBasisGradAtBasisEPoints),
337  basisCurlAtBasisCurlEPoints_(basisCurlAtBasisCurlEPoints),
338  hgradBasisGradAtTargetEPoints_(hgradBasisGradAtTargetEPoints),
339  basisCurlAtTargetCurlEPoints_(basisCurlAtTargetCurlEPoints),
340  basisEWeights_(basisEWeights), basisCurlEWeights_(basisCurlEWeights),
341  wHgradBasisGradAtBasisEPoints_(wHgradBasisGradAtBasisEPoints),
342  wBasisCurlAtBasisCurlEPoints_(wBasisCurlAtBasisCurlEPoints),
343  targetEWeights_(targetEWeights), targetCurlEWeights_(targetCurlEWeights),
344  wHgradBasisGradAtTargetEPoints_(wHgradBasisGradAtTargetEPoints),
345  wBasisCurlAtTargetCurlEPoints_(wBasisCurlAtTargetCurlEPoints),
346  computedDofs_(computedDofs), tagToOrdinal_(tagToOrdinal), hGradTagToOrdinal_(hGradTagToOrdinal),
347  numCellDofs_(numCellDofs), hgradCardinality_(hgradCardinality),
348  offsetBasis_(offsetBasis), offsetBasisCurl_(offsetBasisCurl), offsetTargetCurl_(offsetTargetCurl),
349  numEdgeFaceDofs_(numEdgeFaceDofs), dim_(dim), derDim_(derDim) {}
350 
351  void
352  KOKKOS_INLINE_FUNCTION
353  operator()(const ordinal_type ic) const {
354 
355  ordinal_type numBasisPoints = basisEWeights_.extent(0);
356  ordinal_type numBasisCurlPoints = basisCurlEWeights_.extent(0);
357  ordinal_type numTargetPoints = targetEWeights_.extent(0);
358  ordinal_type numTargetCurlPoints = targetCurlEWeights_.extent(0);
359  for(ordinal_type j=0; j <hgradCardinality_; ++j) {
360  ordinal_type idof = hGradTagToOrdinal_(dim_, 0, j);
361  for(ordinal_type d=0; d <dim_; ++d) {
362  for(ordinal_type iq=0; iq <numBasisPoints; ++iq)
363  wHgradBasisGradAtBasisEPoints_(ic,j,iq,d) = hgradBasisGradAtBasisEPoints_(idof,iq,d)*basisEWeights_(iq);
364  for(ordinal_type iq=0; iq <numTargetPoints; ++iq)
365  wHgradBasisGradAtTargetEPoints_(ic,j,iq,d) = hgradBasisGradAtTargetEPoints_(idof,iq,d)*targetEWeights_(iq);
366  }
367  }
368  for(ordinal_type j=0; j <numCellDofs_; ++j) {
369  ordinal_type idof = tagToOrdinal_(dim_, 0, j);
370  for(ordinal_type d=0; d <dim_; ++d)
371  for(ordinal_type iq=0; iq <numBasisPoints; ++iq)
372  cellBasisAtBasisEPoints_(ic,j,iq,d)=basisAtBasisEPoints_(ic,idof,offsetBasis_+iq,d);
373 
374  for(ordinal_type d=0; d <derDim_; ++d) {
375  for(ordinal_type iq=0; iq <numBasisCurlPoints; ++iq) {
376  cellBasisCurlAtBasisCurlEPoints_(ic,j,iq,d)=basisCurlAtBasisCurlEPoints_(ic,idof,offsetBasisCurl_+iq,d);
377  wBasisCurlAtBasisCurlEPoints_(ic,j,iq,d)=cellBasisCurlAtBasisCurlEPoints_(ic,j,iq,d)*basisCurlEWeights_(iq);
378  }
379  for(ordinal_type iq=0; iq <numTargetCurlPoints; ++iq)
380  wBasisCurlAtTargetCurlEPoints_(ic,j,iq,d) = basisCurlAtTargetCurlEPoints_(ic,idof,offsetTargetCurl_+iq,d)*targetCurlEWeights_(iq);
381  }
382  }
383  for(ordinal_type j=0; j < numEdgeFaceDofs_; ++j) {
384  ordinal_type jdof = computedDofs_(j);
385  for(ordinal_type d=0; d <derDim_; ++d)
386  for(ordinal_type iq=0; iq <numBasisCurlPoints; ++iq)
387  negPartialProjCurl_(ic,iq,d) -= basisCoeffs_(ic,jdof)*basisCurlAtBasisCurlEPoints_(ic,jdof,offsetBasisCurl_+iq,d);
388  for(ordinal_type d=0; d <dim_; ++d)
389  for(ordinal_type iq=0; iq <numBasisPoints; ++iq)
390  negPartialProj_(ic,iq,d) -= basisCoeffs_(ic,jdof)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
391  }
392  }
393 };
394 
395 
396 template<typename SpT>
397 template<typename BasisType,
398 typename ortValueType, class ...ortProperties>
399 void
400 ProjectionTools<SpT>::getHCurlEvaluationPoints(typename BasisType::ScalarViewType targetEPoints,
401  typename BasisType::ScalarViewType targetCurlEPoints,
402  const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
403  const BasisType* cellBasis,
405  const EvalPointsType evalPointType) {
406  typedef typename BasisType::scalarType scalarType;
407  typedef Kokkos::DynRankView<scalarType,ortProperties...> ScalarViewType;
408  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
409  const auto cellTopo = cellBasis->getBaseCellTopology();
410  ordinal_type dim = cellTopo.getDimension();
411  ordinal_type numCells = targetEPoints.extent(0);
412  const ordinal_type edgeDim = 1;
413  const ordinal_type faceDim = 2;
414 
415  ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;
416  ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount() : 0;
417 
419  typename CellTools<SpT>::subcellParamViewType subcellParamEdge, subcellParamFace;
420  if(numEdges>0)
421  CellTools<SpT>::getSubcellParametrization(subcellParamEdge, edgeDim, cellTopo);
422  if(numFaces>0)
423  CellTools<SpT>::getSubcellParametrization(subcellParamFace, faceDim, cellTopo);
424 
425  auto refTopologyKey = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTopologyKey());
426 
427  auto evalPointsRange = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getPointsRange(evalPointType));
428  auto curlEPointsRange = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getDerivPointsRange(evalPointType));
429 
430 
431  ScalarViewType workView("workView", numCells, std::max(projStruct->getMaxNumEvalPoints(evalPointType),projStruct->getMaxNumDerivPoints(evalPointType)), dim-1);
432  for(ordinal_type ie=0; ie<numEdges; ++ie) {
433 
434  auto edgePointsRange = evalPointsRange(edgeDim, ie);
435  auto edgeRefPointsRange = range_type(0, range_size(edgePointsRange));
436  auto edgeEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getEvalPoints(edgeDim,ie,evalPointType));
437 
438  Kokkos::parallel_for
439  ("Evaluate Points Edges ",
440  Kokkos::RangePolicy<SpT, int> (0, numCells),
441  KOKKOS_LAMBDA (const size_t ic) {
442 
443  ordinal_type eOrt[12];
444  orts(ic).getEdgeOrientation(eOrt, numEdges);
445  ordinal_type ort = eOrt[ie];
446 
447  auto orientedEdgeEPoints = Kokkos::subview(workView, ic, edgeRefPointsRange, range_type(0,edgeDim));
448  Impl::OrientationTools::mapToModifiedReference(orientedEdgeEPoints,edgeEPoints,refTopologyKey(edgeDim, ie),ort);
449  CellTools<SpT>::mapToReferenceSubcell(Kokkos::subview(targetEPoints,ic,edgePointsRange,Kokkos::ALL()), orientedEdgeEPoints, subcellParamEdge, edgeDim, ie, dim);
450  });
451  }
452 
453 
454  for(ordinal_type iface=0; iface<numFaces; ++iface) {
455  auto faceCurlPointsRange = curlEPointsRange(faceDim, iface);
456  auto faceRefCurlPointsRange = range_type(0, range_size(faceCurlPointsRange));
457  auto faceCurlEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getDerivEvalPoints(faceDim,iface,evalPointType));
458 
459  auto facePointsRange = evalPointsRange(faceDim, iface);
460  auto faceRefPointsRange = range_type(0, range_size(facePointsRange));
461  auto faceEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getEvalPoints(faceDim,iface,evalPointType));
462 
463  Kokkos::parallel_for
464  ("Evaluate Points Faces ",
465  Kokkos::RangePolicy<SpT, int> (0, numCells),
466  KOKKOS_LAMBDA (const size_t ic) {
467 
468  ordinal_type fOrt[6];
469  orts(ic).getFaceOrientation(fOrt, numFaces);
470  ordinal_type ort = fOrt[iface];
471 
472 
473  auto orientedFaceEPoints = Kokkos::subview(workView, ic, faceRefPointsRange, Kokkos::ALL());
474  Impl::OrientationTools::mapToModifiedReference(orientedFaceEPoints,faceEPoints,refTopologyKey(faceDim, iface),ort);
475  CellTools<SpT>::mapToReferenceSubcell(Kokkos::subview(targetEPoints, ic, facePointsRange, Kokkos::ALL()), orientedFaceEPoints , subcellParamFace, faceDim, iface, dim);
476 
477  auto orientedFaceCurlEPoints = Kokkos::subview(workView, ic, faceRefCurlPointsRange, Kokkos::ALL());
478  Impl::OrientationTools::mapToModifiedReference(orientedFaceCurlEPoints,faceCurlEPoints,refTopologyKey(faceDim, iface),ort);
479  CellTools<SpT>::mapToReferenceSubcell(Kokkos::subview(targetCurlEPoints, ic, faceCurlPointsRange, Kokkos::ALL()), orientedFaceCurlEPoints, subcellParamFace, faceDim, iface, dim);
480  });
481  }
482 
483 
484  if(cellBasis->getDofCount(dim,0)>0) {
485  auto cellEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getEvalPoints(dim,0,evalPointType));
486  RealSpaceTools<SpT>::clone(Kokkos::subview(targetEPoints, Kokkos::ALL(), evalPointsRange(dim, 0), Kokkos::ALL()), cellEPoints);
487 
488  auto cellCurlEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getDerivEvalPoints(dim,0,evalPointType));
489  RealSpaceTools<SpT>::clone(Kokkos::subview(targetCurlEPoints, Kokkos::ALL(), curlEPointsRange(dim, 0), Kokkos::ALL()), cellCurlEPoints);
490  }
491 }
492 
493 
494 template<typename SpT>
495 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
496 typename funValsValueType, class ...funValsProperties,
497 typename BasisType,
498 typename ortValueType,class ...ortProperties>
499 void
500 ProjectionTools<SpT>::getHCurlBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
501  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtTargetEPoints,
502  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetCurlAtTargetCurlEPoints,
503  const typename BasisType::ScalarViewType targetEPoints,
504  const typename BasisType::ScalarViewType targetCurlEPoints,
505  const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
506  const BasisType* cellBasis,
507  ProjectionStruct<SpT, typename BasisType::scalarType> * projStruct){
508 
509  typedef typename BasisType::scalarType scalarType;
510  typedef Kokkos::DynRankView<scalarType,SpT> ScalarViewType;
511  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
512  const auto cellTopo = cellBasis->getBaseCellTopology();
513  ordinal_type dim = cellTopo.getDimension();
514  ordinal_type numTotalTargetEPoints(targetAtTargetEPoints.extent(1)),
515  numTotalTargetCurlEPoints(targetCurlAtTargetCurlEPoints.extent(1));
516  ordinal_type basisCardinality = cellBasis->getCardinality();
517  ordinal_type numCells = targetAtTargetEPoints.extent(0);
518  const ordinal_type edgeDim = 1;
519  const ordinal_type faceDim = 2;
520  const ordinal_type derDim = dim == 3 ? dim : 1;
521 
522  const Kokkos::RangePolicy<SpT> policy(0, numCells);
523 
524  const std::string& name = cellBasis->getName();
525 
526  ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;
527  ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount() : 0;
528 
529  ScalarViewType refEdgesTangent("refEdgesTangent", numEdges, dim);
530  ScalarViewType refFacesTangents("refFaceTangents", numFaces, dim, 2);
531  ScalarViewType refFacesNormal("refFaceNormal", numFaces, dim);
532 
533  ordinal_type numEdgeDofs(0);
534  for(ordinal_type ie=0; ie<numEdges; ++ie)
535  numEdgeDofs += cellBasis->getDofCount(edgeDim,ie);
536 
537  ordinal_type numFaceDofs(0);
538  for(ordinal_type iface=0; iface<numFaces; ++iface)
539  numFaceDofs += cellBasis->getDofCount(faceDim,iface);
540 
541  auto tagToOrdinal = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(), cellBasis->getAllDofOrdinal());
542 
543  Kokkos::View<ordinal_type*, SpT> computedDofs("computedDofs",numEdgeDofs+numFaceDofs);
544 
545  auto targetEPointsRange = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetPointsRange());
546  auto targetCurlEPointsRange = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetDerivPointsRange());
547 
548  auto basisEPointsRange = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisPointsRange());
549  auto basisCurlEPointsRange = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisDerivPointsRange());
550 
551  auto refTopologyKey = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTopologyKey());
552 
553  ordinal_type numTotalBasisEPoints = projStruct->getNumBasisEvalPoints(), numTotalBasisCurlEPoints = projStruct->getNumBasisDerivEvalPoints();
554 
555  ScalarViewType basisEPoints("basisEPoints",numCells,numTotalBasisEPoints, dim);
556  ScalarViewType basisCurlEPoints("basisCurlEPoints",numCells,numTotalBasisCurlEPoints, dim);
557  getHCurlEvaluationPoints(basisEPoints, basisCurlEPoints, orts, cellBasis, projStruct, EvalPointsType::BASIS);
558 
559  ScalarViewType basisAtBasisEPoints("basisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints, dim);
560  ScalarViewType basisAtTargetEPoints("basisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints, dim);
561  {
562  ScalarViewType nonOrientedBasisAtBasisEPoints("nonOrientedBasisAtEPoints",numCells,basisCardinality, numTotalBasisEPoints, dim);
563  ScalarViewType nonOrientedBasisAtTargetEPoints("nonOrientedBasisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints, dim);
564  for(ordinal_type ic=0; ic<numCells; ++ic) {
565  cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtTargetEPoints,ic,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(targetEPoints, ic, Kokkos::ALL(), Kokkos::ALL()));
566  cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtBasisEPoints,ic,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(basisEPoints, ic, Kokkos::ALL(), Kokkos::ALL()));
567  }
568 
569  OrientationTools<SpT>::modifyBasisByOrientation(basisAtBasisEPoints, nonOrientedBasisAtBasisEPoints, orts, cellBasis);
570  OrientationTools<SpT>::modifyBasisByOrientation(basisAtTargetEPoints, nonOrientedBasisAtTargetEPoints, orts, cellBasis);
571  }
572 
573  ScalarViewType basisCurlAtBasisCurlEPoints;
574  ScalarViewType basisCurlAtTargetCurlEPoints;
575  if(numTotalBasisCurlEPoints>0) {
576  ScalarViewType nonOrientedBasisCurlAtTargetCurlEPoints, nonOrientedBasisCurlAtBasisCurlEPoints;
577  if (dim == 3) {
578  basisCurlAtBasisCurlEPoints = ScalarViewType ("basisCurlAtBasisCurlEPoints",numCells,basisCardinality, numTotalBasisCurlEPoints, dim);
579  nonOrientedBasisCurlAtBasisCurlEPoints = ScalarViewType ("nonOrientedBasisCurlAtBasisCurlEPoints",numCells,basisCardinality, numTotalBasisCurlEPoints, dim);
580  basisCurlAtTargetCurlEPoints = ScalarViewType("basisCurlAtTargetCurlEPoints",numCells,basisCardinality, numTotalTargetCurlEPoints, dim);
581  nonOrientedBasisCurlAtTargetCurlEPoints = ScalarViewType("nonOrientedBasisCurlAtTargetCurlEPoints",numCells,basisCardinality, numTotalTargetCurlEPoints, dim);
582  } else {
583  basisCurlAtBasisCurlEPoints = ScalarViewType ("basisCurlAtBasisCurlEPoints",numCells,basisCardinality, numTotalBasisCurlEPoints);
584  nonOrientedBasisCurlAtBasisCurlEPoints = ScalarViewType ("nonOrientedBasisCurlAtBasisCurlEPoints",numCells,basisCardinality, numTotalBasisCurlEPoints);
585  basisCurlAtTargetCurlEPoints = ScalarViewType("basisCurlAtTargetCurlEPoints",numCells,basisCardinality, numTotalTargetCurlEPoints);
586  nonOrientedBasisCurlAtTargetCurlEPoints = ScalarViewType("nonOrientedBasisCurlAtTargetCurlEPoints",numCells,basisCardinality, numTotalTargetCurlEPoints);
587  }
588  for(ordinal_type ic=0; ic<numCells; ++ic) {
589  cellBasis->getValues(Kokkos::subview(nonOrientedBasisCurlAtBasisCurlEPoints,ic,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(basisCurlEPoints, ic, Kokkos::ALL(), Kokkos::ALL()),OPERATOR_CURL);
590  cellBasis->getValues(Kokkos::subview(nonOrientedBasisCurlAtTargetCurlEPoints,ic,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(targetCurlEPoints, ic, Kokkos::ALL(), Kokkos::ALL()),OPERATOR_CURL);
591  }
592  OrientationTools<SpT>::modifyBasisByOrientation(basisCurlAtBasisCurlEPoints, nonOrientedBasisCurlAtBasisCurlEPoints, orts, cellBasis);
593  OrientationTools<SpT>::modifyBasisByOrientation(basisCurlAtTargetCurlEPoints, nonOrientedBasisCurlAtTargetCurlEPoints, orts, cellBasis);
594  }
595 
596  ordinal_type computedDofsCount = 0;
597  for(ordinal_type ie=0; ie<numEdges; ++ie) {
598 
599  ordinal_type edgeCardinality = cellBasis->getDofCount(edgeDim,ie);
600  ordinal_type numBasisEPoints = range_size(basisEPointsRange(edgeDim, ie));
601  ordinal_type numTargetEPoints = range_size(targetEPointsRange(edgeDim, ie));
602 
603  {
604  auto refEdgeTan = Kokkos::subview(refEdgesTangent, ie, Kokkos::ALL());
605  auto refEdgeTanHost = Kokkos::create_mirror_view(refEdgeTan);
606  CellTools<SpT>::getReferenceEdgeTangent(refEdgeTanHost, ie, cellTopo);
607  Kokkos::deep_copy(refEdgeTan,refEdgeTanHost);
608  }
609 
610  ScalarViewType basisTanAtBasisEPoints("basisTanAtBasisEPoints",numCells,edgeCardinality, numBasisEPoints);
611  ScalarViewType basisTanAtTargetEPoints("basisTanAtTargetEPoints",numCells,edgeCardinality, numTargetEPoints);
612  ScalarViewType weightedTanBasisAtBasisEPoints("weightedTanBasisAtBasisEPoints",numCells,edgeCardinality, numBasisEPoints);
613  ScalarViewType weightedTanBasisAtTargetEPoints("weightedTanBasisAtTargetEPoints",numCells,edgeCardinality, numTargetEPoints);
614  ScalarViewType targetTanAtTargetEPoints("normalTargetAtTargetEPoints",numCells, numTargetEPoints);
615 
616  auto targetEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetEvalWeights(edgeDim,ie));
617  auto basisEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisEvalWeights(edgeDim,ie));
618 
619  //Note: we are not considering the jacobian of the orientation map since it is simply a scalar term for the integrals and it does not affect the projection
620  ordinal_type offsetBasis = basisEPointsRange(edgeDim, ie).first;
621  ordinal_type offsetTarget = targetEPointsRange(edgeDim, ie).first;
622 
623  typedef ComputeBasisCoeffsOnEdges_HCurl<ScalarViewType, decltype(basisEWeights), decltype(tagToOrdinal), decltype(targetAtTargetEPoints)> functorTypeEdge;
624  Kokkos::parallel_for(policy, functorTypeEdge(basisTanAtBasisEPoints,basisAtBasisEPoints,basisEWeights,
625  weightedTanBasisAtBasisEPoints, targetEWeights,
626  basisAtTargetEPoints, weightedTanBasisAtTargetEPoints, tagToOrdinal,
627  targetAtTargetEPoints, targetTanAtTargetEPoints,
628  refEdgesTangent, edgeCardinality, offsetBasis,
629  offsetTarget, edgeDim,
630  dim, ie));
631 
632  ScalarViewType edgeMassMat_("edgeMassMat_", numCells, edgeCardinality+1, edgeCardinality+1),
633  edgeRhsMat_("rhsMat_", numCells, edgeCardinality+1);
634 
635  ScalarViewType eWeights_("eWeights_", numCells, 1, basisEWeights.extent(0)), targetEWeights_("targetEWeights", numCells, 1, targetEWeights.extent(0));
636  RealSpaceTools<SpT>::clone(eWeights_, basisEWeights);
637  RealSpaceTools<SpT>::clone(targetEWeights_, targetEWeights);
638 
639  range_type range_H(0, edgeCardinality);
640  range_type range_B(edgeCardinality, edgeCardinality+1);
641  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(edgeMassMat_,Kokkos::ALL(),range_H,range_H), basisTanAtBasisEPoints, weightedTanBasisAtBasisEPoints);
642  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(edgeMassMat_,Kokkos::ALL(),range_H,range_B), basisTanAtBasisEPoints, eWeights_);
643  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(edgeRhsMat_,Kokkos::ALL(),range_H), targetTanAtTargetEPoints, weightedTanBasisAtTargetEPoints);
644  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(edgeRhsMat_,Kokkos::ALL(),range_B), targetTanAtTargetEPoints, targetEWeights_);
645 
646  typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, SpT> WorkArrayViewType;
647  ScalarViewType t_("t",numCells, edgeCardinality+1);
648  WorkArrayViewType w_("w",numCells, edgeCardinality+1);
649 
650  auto edgeDofs = Kokkos::subview(tagToOrdinal, edgeDim, ie, Kokkos::ALL());
651  ElemSystem edgeSystem("edgeSystem", false);
652  edgeSystem.solve(basisCoeffs, edgeMassMat_, edgeRhsMat_, t_, w_, edgeDofs, edgeCardinality, 1);
653 
654  for(ordinal_type i=0; i<edgeCardinality; ++i)
655  computedDofs(computedDofsCount++) = cellBasis->getDofOrdinal(edgeDim, ie, i);
656 
657  }
658 
659  Basis<SpT,scalarType,scalarType> *hgradBasis = NULL;
660  for(ordinal_type iface=0; iface<numFaces; ++iface) {
661 
662  if(cellTopo.getKey() == shards::getCellTopologyData<shards::Hexahedron<8> >()->key)
663  hgradBasis = new Basis_HGRAD_QUAD_Cn_FEM<SpT,scalarType,scalarType>(cellBasis->getDegree(),POINTTYPE_WARPBLEND);
664  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Tetrahedron<4> >()->key)
665  hgradBasis = new Basis_HGRAD_TRI_Cn_FEM<SpT,scalarType,scalarType>(cellBasis->getDegree(),POINTTYPE_WARPBLEND);
666  else {
667  std::stringstream ss;
668  ss << ">>> ERROR (Intrepid2::ProjectionTools::getHCurlBasisCoeffs): "
669  << "Method not implemented for basis " << name;
670  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
671  }
672 
673  {
674  auto refFaceTanU = Kokkos::subview(refFacesTangents, iface, Kokkos::ALL, 0);
675  auto refFaceTanV = Kokkos::subview(refFacesTangents, iface, Kokkos::ALL,1);
676  auto refFaceTanUHost = Kokkos::create_mirror_view(refFaceTanU);
677  auto refFaceTanVHost = Kokkos::create_mirror_view(refFaceTanV);
678  CellTools<SpT>::getReferenceFaceTangents(refFaceTanUHost, refFaceTanVHost, iface, cellTopo);
679  Kokkos::deep_copy(refFaceTanU, refFaceTanUHost);
680  Kokkos::deep_copy(refFaceTanV, refFaceTanVHost);
681  auto refFaceNormal = Kokkos::subview(refFacesNormal,iface,Kokkos::ALL());
682  auto refFaceNormalHost = Kokkos::create_mirror_view(refFaceNormal);
683  CellTools<SpT>::getReferenceFaceNormal(refFaceNormalHost, iface, cellTopo);
684  Kokkos::deep_copy(refFaceNormal, refFaceNormalHost);
685  }
686 
687  ordinal_type numTargetEPoints = range_size(targetEPointsRange(faceDim, iface));
688  ordinal_type numTargetCurlEPoints = range_size(targetCurlEPointsRange(faceDim, iface));
689  ordinal_type numBasisEPoints = range_size(basisEPointsRange(faceDim, iface));
690  ordinal_type numBasisCurlEPoints = range_size(basisCurlEPointsRange(faceDim, iface));
691 
692  ordinal_type numFaceDofs = cellBasis->getDofCount(faceDim,iface);
693 
694  ScalarViewType hgradBasisGradAtBasisEPoints("hgradBasisGradAtBasisEPoints",hgradBasis->getCardinality(), numBasisEPoints, faceDim);
695  ScalarViewType hgradBasisGradAtTargetEPoints("hgradBasisGradAtTargetEPoints",hgradBasis->getCardinality(), numTargetEPoints, faceDim);
696 
697  ordinal_type hgradCardinality = hgradBasis->getDofCount(faceDim,0);
698 
699  auto refBasisEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisEvalPoints(faceDim, iface));
700  auto refTargetEPoints = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetEvalPoints(faceDim, iface));
701  hgradBasis->getValues(hgradBasisGradAtBasisEPoints, refBasisEPoints, OPERATOR_GRAD);
702  hgradBasis->getValues(hgradBasisGradAtTargetEPoints, refTargetEPoints, OPERATOR_GRAD);
703 
704  ScalarViewType basisTanAtBasisEPoints("basisTanAtBasisEPoints",numCells,numFaceDofs, numBasisEPoints,dim-1);
705  ScalarViewType basisTanAtTargetEPoints("basisTanAtTargetEPoints",numCells,numFaceDofs, numTargetEPoints,dim-1);
706  ScalarViewType basisCurlNormalAtBasisCurlEPoints("normaBasisCurlAtBasisEPoints",numCells,numFaceDofs, numBasisCurlEPoints);
707  ScalarViewType wNormalBasisCurlAtBasisCurlEPoints("weightedNormalBasisCurlAtBasisEPoints",numCells,numFaceDofs, numBasisCurlEPoints);
708 
709  ScalarViewType targetTanAtTargetEPoints("targetTanAtTargetEPoints",numCells, numTargetEPoints, dim-1);
710  ScalarViewType normalTargetCurlAtTargetEPoints("normalTargetCurlAtTargetEPoints",numCells, numTargetCurlEPoints);
711  ScalarViewType wNormalBasisCurlBasisAtTargetCurlEPoints("weightedNormalBasisCurlAtTargetCurlEPoints",numCells,numFaceDofs, numTargetCurlEPoints);
712 
713  ScalarViewType wHgradBasisGradAtBasisEPoints("wHgradBasisGradAtBasisEPoints",numCells, hgradCardinality, numBasisEPoints, faceDim);
714  ScalarViewType wHgradBasisGradAtTargetEPoints("wHgradBasisGradAtTargetEPoints",numCells, hgradCardinality, numTargetEPoints, faceDim);
715 
716  ScalarViewType negPartialProjCurlNormal("mNormalComputedProjection", numCells,numBasisEPoints);
717  ScalarViewType negPartialProjTan("negPartialProjTan", numCells,numBasisEPoints,dim-1);
718 
719 
720  ordinal_type offsetBasis = basisEPointsRange(faceDim, iface).first;
721  ordinal_type offsetBasisCurl = basisCurlEPointsRange(faceDim, iface).first;
722  ordinal_type offsetTarget = targetEPointsRange(faceDim, iface).first;
723  ordinal_type offsetTargetCurl = targetCurlEPointsRange(faceDim, iface).first;
724 
725  auto hGradTagToOrdinal = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(), hgradBasis->getAllDofOrdinal());
726 
727  auto basisEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisEvalWeights(faceDim,iface));
728  auto targetEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetEvalWeights(faceDim,iface));
729  auto targetCurlEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetDerivEvalWeights(faceDim,iface));
730  auto basisCurlEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisDerivEvalWeights(faceDim,iface));
731  typedef ComputeBasisCoeffsOnFaces_HCurl<decltype(basisCoeffs), decltype(orts), ScalarViewType, decltype(targetAtTargetEPoints), decltype(basisEWeights),
732  decltype(tagToOrdinal), decltype(refTopologyKey), decltype(computedDofs)> functorTypeFaces;
733  Kokkos::parallel_for(policy, functorTypeFaces(basisCoeffs,
734  orts, negPartialProjTan, negPartialProjCurlNormal,
735  hgradBasisGradAtBasisEPoints, wHgradBasisGradAtBasisEPoints,
736  basisCurlAtBasisCurlEPoints, basisCurlNormalAtBasisCurlEPoints,
737  basisAtBasisEPoints,
738  normalTargetCurlAtTargetEPoints, basisTanAtBasisEPoints,
739  hgradBasisGradAtTargetEPoints, wHgradBasisGradAtTargetEPoints,
740  wNormalBasisCurlAtBasisCurlEPoints, basisCurlAtTargetCurlEPoints,
741  wNormalBasisCurlBasisAtTargetCurlEPoints, targetAtTargetEPoints,
742  targetTanAtTargetEPoints, targetCurlAtTargetCurlEPoints,
743  basisEWeights, targetEWeights,
744  basisCurlEWeights, targetCurlEWeights, tagToOrdinal,
745  hGradTagToOrdinal, refTopologyKey,
746  refFacesNormal, refFacesTangents,
747  computedDofs, offsetBasis,
748  offsetBasisCurl, offsetTarget,
749  offsetTargetCurl, iface,
750  hgradCardinality, numFaces,
751  numFaceDofs, numEdgeDofs,
752  faceDim, dim));
753 
754 
755  ScalarViewType faceMassMat_("faceMassMat_", numCells, numFaceDofs+hgradCardinality, numFaceDofs+hgradCardinality),
756  faceRhsMat_("rhsMat_", numCells, numFaceDofs+hgradCardinality);
757  range_type range_H(0, numFaceDofs);
758  range_type range_B(numFaceDofs, numFaceDofs+hgradCardinality);
759  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(faceMassMat_,Kokkos::ALL(),range_H,range_H), basisCurlNormalAtBasisCurlEPoints, wNormalBasisCurlAtBasisCurlEPoints);
760  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(faceMassMat_,Kokkos::ALL(),range_H,range_B), basisTanAtBasisEPoints, wHgradBasisGradAtBasisEPoints);
761 
762  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(faceRhsMat_,Kokkos::ALL(),range_H), normalTargetCurlAtTargetEPoints, wNormalBasisCurlBasisAtTargetCurlEPoints);
763  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(faceRhsMat_,Kokkos::ALL(),range_H), negPartialProjCurlNormal, wNormalBasisCurlAtBasisCurlEPoints,true);
764 
765  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(faceRhsMat_,Kokkos::ALL(),range_B), targetTanAtTargetEPoints, wHgradBasisGradAtTargetEPoints);
766  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(faceRhsMat_,Kokkos::ALL(),range_B), negPartialProjTan, wHgradBasisGradAtBasisEPoints,true);
767 
768 
769  typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, SpT> WorkArrayViewType;
770  ScalarViewType t_("t",numCells, numFaceDofs+hgradCardinality);
771  WorkArrayViewType w_("w",numCells, numFaceDofs+hgradCardinality);
772 
773  auto faceDofs = Kokkos::subview(tagToOrdinal, faceDim, iface, Kokkos::ALL());
774  ElemSystem faceSystem( "faceSystem", false);
775  faceSystem.solve(basisCoeffs, faceMassMat_, faceRhsMat_, t_, w_, faceDofs, numFaceDofs, hgradCardinality);
776 
777  for(ordinal_type i=0; i<numFaceDofs; ++i)
778  computedDofs(computedDofsCount++) = cellBasis->getDofOrdinal(faceDim, iface, i);
779 
780  delete hgradBasis;
781  }
782 
783  ordinal_type numCellDofs = cellBasis->getDofCount(dim,0);
784  if(numCellDofs>0) {
785  if(cellTopo.getKey() == shards::getCellTopologyData<shards::Hexahedron<8> >()->key)
786  hgradBasis = new Basis_HGRAD_HEX_Cn_FEM<SpT,scalarType,scalarType>(cellBasis->getDegree());
787  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Tetrahedron<4> >()->key)
788  hgradBasis = new Basis_HGRAD_TET_Cn_FEM<SpT,scalarType,scalarType>(cellBasis->getDegree(),POINTTYPE_WARPBLEND);
789  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Triangle<3> >()->key)
790  hgradBasis = new Basis_HGRAD_TRI_Cn_FEM<SpT,scalarType,scalarType>(cellBasis->getDegree(),POINTTYPE_WARPBLEND);
791  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Quadrilateral<4> >()->key)
792  hgradBasis = new Basis_HGRAD_QUAD_Cn_FEM<SpT,scalarType,scalarType>(cellBasis->getDegree(),POINTTYPE_WARPBLEND);
793  else {
794  std::stringstream ss;
795  ss << ">>> ERROR (Intrepid2::ProjectionTools::getHCurlBasisCoeffs): "
796  << "Method not implemented for basis " << name;
797  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
798  }
799 
800  range_type cellPointsRange = targetEPointsRange(dim, 0);
801  range_type cellCurlPointsRange = targetCurlEPointsRange(dim, 0);
802 
803  ordinal_type numTargetCurlEPoints = range_size(targetCurlEPointsRange(dim,0));
804  ordinal_type numBasisCurlEPoints = range_size(basisCurlEPointsRange(dim,0));
805  ordinal_type numBasisEPoints = range_size(basisEPointsRange(dim,0));
806  ordinal_type numTargetEPoints = range_size(targetEPointsRange(dim,0));
807 
808  ScalarViewType hgradBasisGradAtBasisEPoints("hgradBasisGradAtBasisEPoints",hgradBasis->getCardinality(), numBasisEPoints, dim);
809  ScalarViewType hgradBasisGradAtTargetEPoints("hgradBasisGradAtTargetEPoints",hgradBasis->getCardinality(), numTargetEPoints, dim);
810 
811  ordinal_type hgradCardinality = hgradBasis->getDofCount(dim,0);
812  ScalarViewType wHgradBasisGradAtTargetEPoints("wHgradBasisGradAtTargetEPoints",numCells, hgradCardinality, numTargetEPoints, dim);
813  ScalarViewType wHgradBasisGradAtBasisEPoints("wHgradBasisGradAtBasisEPoints",numCells, hgradCardinality, numBasisEPoints, dim);
814 
815  hgradBasis->getValues(hgradBasisGradAtBasisEPoints,Kokkos::subview(basisEPoints, 0, basisEPointsRange(dim, 0), Kokkos::ALL()), OPERATOR_GRAD);
816  hgradBasis->getValues(hgradBasisGradAtTargetEPoints,Kokkos::subview(targetEPoints, 0, targetEPointsRange(dim, 0), Kokkos::ALL()),OPERATOR_GRAD);
817 
818  ScalarViewType cellBasisAtBasisEPoints("basisCellAtEPoints",numCells,numCellDofs, numBasisEPoints, dim);
819  ScalarViewType cellBasisCurlAtCurlEPoints("cellBasisCurlAtCurlEPoints",numCells,numCellDofs, numBasisCurlEPoints, derDim);
820  ScalarViewType negPartialProjCurl("negPartialProjCurl", numCells, numBasisEPoints, derDim);
821  ScalarViewType negPartialProj("negPartialProj", numCells, numBasisEPoints, dim);
822  ScalarViewType wBasisCurlAtCurlEPoints("weightedBasisCurlAtBasisEPoints",numCells,numCellDofs, numBasisCurlEPoints,derDim);
823  ScalarViewType wBasisCurlBasisAtTargetCurlEPoints("weightedBasisCurlAtTargetCurlEPoints",numCells,numCellDofs, numTargetCurlEPoints,derDim);
824 
825  auto targetEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetEvalWeights(dim,0));
826  auto basisEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisEvalWeights(dim,0));
827  auto targetCurlEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getTargetDerivEvalWeights(dim,0));
828  auto basisCurlEWeights = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(),projStruct->getBasisDerivEvalWeights(dim,0));
829  ordinal_type offsetBasis = basisEPointsRange(dim, 0).first;
830  ordinal_type offsetBasisCurl = basisCurlEPointsRange(dim, 0).first;
831  ordinal_type offsetTargetCurl = targetCurlEPointsRange(dim, 0).first;
832 
833 
834  auto hGradTagToOrdinal = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(), hgradBasis->getAllDofOrdinal());
835 
836  typedef ComputeBasisCoeffsOnCell_HCurl<decltype(basisCoeffs), ScalarViewType, decltype(basisEWeights),
837  decltype(computedDofs), decltype(tagToOrdinal)> functorTypeCell;
838  Kokkos::parallel_for(policy, functorTypeCell(basisCoeffs, negPartialProj, negPartialProjCurl,
839  cellBasisAtBasisEPoints, cellBasisCurlAtCurlEPoints,
840  basisAtBasisEPoints, hgradBasisGradAtBasisEPoints, basisCurlAtBasisCurlEPoints,
841  hgradBasisGradAtTargetEPoints, basisCurlAtTargetCurlEPoints,
842  basisEWeights, basisCurlEWeights, wHgradBasisGradAtBasisEPoints,
843  wBasisCurlAtCurlEPoints, targetEWeights, targetCurlEWeights,
844  wHgradBasisGradAtTargetEPoints,
845  wBasisCurlBasisAtTargetCurlEPoints, computedDofs,
846  tagToOrdinal, hGradTagToOrdinal,
847  numCellDofs, hgradCardinality,
848  offsetBasis, offsetBasisCurl, offsetTargetCurl,
849  numEdgeDofs+numFaceDofs, dim, derDim));
850 
851  ScalarViewType cellMassMat_("cellMassMat_", numCells, numCellDofs+hgradCardinality, numCellDofs+hgradCardinality),
852  cellRhsMat_("rhsMat_", numCells, numCellDofs+hgradCardinality);
853 
854  range_type range_H(0, numCellDofs);
855  range_type range_B(numCellDofs, numCellDofs+hgradCardinality);
856  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellMassMat_,Kokkos::ALL(),range_H,range_H), cellBasisCurlAtCurlEPoints, wBasisCurlAtCurlEPoints);
857  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellMassMat_,Kokkos::ALL(),range_H,range_B), cellBasisAtBasisEPoints, wHgradBasisGradAtBasisEPoints);
858  if(dim==3)
859  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellRhsMat_,Kokkos::ALL(),range_H), Kokkos::subview(targetCurlAtTargetCurlEPoints,Kokkos::ALL(),cellCurlPointsRange,Kokkos::ALL()), wBasisCurlBasisAtTargetCurlEPoints);
860  else
861  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellRhsMat_,Kokkos::ALL(),range_H), Kokkos::subview(targetCurlAtTargetCurlEPoints,Kokkos::ALL(),cellCurlPointsRange), Kokkos::subview(wBasisCurlBasisAtTargetCurlEPoints,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),0));
862  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellRhsMat_,Kokkos::ALL(),range_H), negPartialProjCurl, wBasisCurlAtCurlEPoints, true);
863 
864  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellRhsMat_,Kokkos::ALL(),range_B), Kokkos::subview(targetAtTargetEPoints,Kokkos::ALL(),cellPointsRange,Kokkos::ALL()), wHgradBasisGradAtTargetEPoints);
865  FunctionSpaceTools<SpT >::integrate(Kokkos::subview(cellRhsMat_,Kokkos::ALL(),range_B), negPartialProj, wHgradBasisGradAtBasisEPoints, true);
866 
867  typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, SpT> WorkArrayViewType;
868  ScalarViewType t_("t",numCells, numCellDofs+hgradCardinality);
869  WorkArrayViewType w_("w",numCells, numCellDofs+hgradCardinality);
870 
871  auto cellDofs = Kokkos::subview(tagToOrdinal, dim, 0, Kokkos::ALL());
872  ElemSystem cellSystem( "cellSystem", true);
873  cellSystem.solve(basisCoeffs, cellMassMat_, cellRhsMat_, t_, w_, cellDofs, numCellDofs, hgradCardinality);
874 
875  delete hgradBasis;
876  }
877 }
878 }
879 }
880 
881 #endif
882 
static void getReferenceEdgeTangent(Kokkos::DynRankView< refEdgeTangentValueType, refEdgeTangentProperties...> refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
Header file for the abstract base class Intrepid2::DefaultCubatureFactory.
const range_tag getPointsRange(const EvalPointsType type) const
Returns the range tag of the basis/target evaluation points in subcells.
const range_tag getDerivPointsRange(const EvalPointsType type) const
Returns the range tag of the basis/target derivative evaluation points on subcells.
static void getJacobianOfOrientationMap(JacobianViewType jacobian, const shards::CellTopology cellTopo, const ordinal_type cellOrt)
Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation.
static void modifyBasisByOrientation(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input, const Kokkos::DynRankView< ortValueType, ortProperties...> orts, const BasisType *basis)
Modify basis due to orientation.
static void integrate(Kokkos::DynRankView< outputValueValueType, outputValueProperties...> outputValues, const Kokkos::DynRankView< leftValueValueType, leftValueProperties...> leftValues, const Kokkos::DynRankView< rightValueValueType, rightValueProperties...> rightValues, const bool sumInto=false)
Contracts leftValues and rightValues arrays on the point and possibly space dimensions and stores the...
static void mapToReferenceSubcell(Kokkos::DynRankView< refSubcellPointValueType, refSubcellPointProperties...> refSubcellPoints, const Kokkos::DynRankView< paramPointValueType, paramPointProperties...> paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
ordinal_type getMaxNumEvalPoints(const EvalPointsType type) const
Returns the maximum number of evaluation points across all the subcells.
view_type getDerivEvalPoints(const ordinal_type subCellDim, const ordinal_type subCellId, EvalPointsType type) const
Returns the evaluation points for basis/target derivatives on a subcell.
Header file for the Intrepid2::FunctionSpaceTools class.
static void getHCurlEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType curlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for HCurl projection.
static void getHCurlBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetCurlAtCurlEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType curlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< ExecSpaceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HCurl projection of the target function.
ordinal_type getMaxNumDerivPoints(const EvalPointsType type) const
Returns the maximum number of derivative evaluation points across all the subcells.
static void clone(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Clone input array.
Header file for Intrepid2::ArrayTools class providing utilities for array operations.
static void getReferenceFaceTangents(Kokkos::DynRankView< refFaceTanUValueType, refFaceTanUProperties...> refFaceTanU, Kokkos::DynRankView< refFaceTanVValueType, refFaceTanVProperties...> refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
static void mapToModifiedReference(outPointViewType outPoints, const refPointViewType refPoints, const shards::CellTopology cellTopo, const ordinal_type cellOrt=0)
Computes modified parameterization maps of 1- and 2-subcells with orientation.
static void setSubcellParametrization()
Defines orientation-preserving parametrizations of reference edges and faces of cell topologies with ...
An helper class to compute the evaluation points and weights needed for performing projections...
static void getSubcellParametrization(subcellParamViewType &subcellParam, const ordinal_type subcellDim, const shards::CellTopology parentCell)
Returns array with the coefficients of the parametrization maps for the edges or faces of a reference...
const key_tag getTopologyKey() const
Returns the key tag for subcells.
static void getReferenceFaceNormal(Kokkos::DynRankView< refFaceNormalValueType, refFaceNormalProperties...> refFaceNormal, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to faces of 3D reference cell.
view_type getEvalPoints(const ordinal_type subCellDim, const ordinal_type subCellId, EvalPointsType type) const
Returns the basis/target evaluation points on a subcell.