Intrepid2
Intrepid2_LagrangianInterpolationDef.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_LAGRANGIANINTERPOLATIONDEF_HPP__
49 #define __INTREPID2_LAGRANGIANINTERPOLATIONDEF_HPP__
50 
52 #include "Intrepid2_ArrayTools.hpp"
54 
55 
56 namespace Intrepid2 {
57 namespace Experimental {
58 
59 
60 
61 template<typename ScalarViewType,
62 typename ortViewType,
63 typename t2oViewType,
64 typename subcellParamViewType,
65 typename intViewType>
67  typedef typename ScalarViewType::value_type value_type;
68 
69  ScalarViewType dofCoords_, dofCoeffs_;
70  const ortViewType orts_;
71  const t2oViewType tagToOrdinal_;
72  const subcellParamViewType edgeParam_, faceParam_;
73  const intViewType edgeInternalDofOrdinals_, facesInternalDofOrdinals_;
74  const ScalarViewType edgeInternalDofCoords_, edgeDofCoeffs_, ortJacobianEdge_, refEdgesTan_, refEdgesNormal_;
75  const ScalarViewType facesInternalDofCoords_, faceDofCoeffs_, ortJacobianFace_, refFaceTangents_, refFacesNormal_;
76  ScalarViewType edgeWorkView_, faceWorkView_;
77  const ordinal_type cellDim_, numEdges_, numFaces_;
78  const ordinal_type edgeTopoKey_, numEdgeInternalDofs_;
79  const intViewType faceTopoKey_, numFacesInternalDofs_;
80  const value_type edgeScale_, faceScale_;
81  const bool isBasisHCURL_, isBasisHDIV_;
82 
83  computeDofCoordsAndCoeffs( ScalarViewType dofCoords,
84  ScalarViewType dofCoeffs,
85  const ortViewType orts,
86  const t2oViewType tagToOrdinal,
87  const subcellParamViewType edgeParam,
88  const subcellParamViewType faceParam,
89  const intViewType edgeInternalDofOrdinals,
90  const intViewType facesInternalDofOrdinals,
91  const ScalarViewType edgeInternalDofCoords,
92  const ScalarViewType edgeDofCoeffs,
93  const ScalarViewType refEdgesTan,
94  const ScalarViewType refEdgesNormal,
95  const ScalarViewType facesInternalDofCoords,
96  const ScalarViewType faceDofCoeffs,
97  const ScalarViewType refFaceTangents,
98  const ScalarViewType refFacesNormal,
99  const ordinal_type cellDim,
100  const ordinal_type numEdges,
101  const ordinal_type numFaces,
102  const ordinal_type edgeTopoKey,
103  const ordinal_type numEdgeInternalDofs,
104  const intViewType faceTopoKey,
105  const intViewType numFacesInternalDofs,
106  const value_type edgeScale,
107  const value_type faceScale,
108  const bool isBasisHCURL,
109  const bool isBasisHDIV
110  )
111  : dofCoords_(dofCoords),
112  dofCoeffs_(dofCoeffs),
113  orts_(orts),
114  tagToOrdinal_(tagToOrdinal),
115  edgeParam_(edgeParam),
116  faceParam_(faceParam),
117  edgeInternalDofOrdinals_(edgeInternalDofOrdinals),
118  facesInternalDofOrdinals_(facesInternalDofOrdinals),
119  edgeInternalDofCoords_(edgeInternalDofCoords),
120  edgeDofCoeffs_(edgeDofCoeffs),
121  refEdgesTan_(refEdgesTan),
122  refEdgesNormal_(refEdgesNormal),
123  facesInternalDofCoords_(facesInternalDofCoords),
124  faceDofCoeffs_(faceDofCoeffs),
125  refFaceTangents_(refFaceTangents),
126  refFacesNormal_(refFacesNormal),
127  cellDim_(cellDim),
128  numEdges_(numEdges),
129  numFaces_(numFaces),
130  edgeTopoKey_(edgeTopoKey),
131  numEdgeInternalDofs_(numEdgeInternalDofs),
132  faceTopoKey_(faceTopoKey),
133  numFacesInternalDofs_(numFacesInternalDofs),
134  edgeScale_(edgeScale),
135  faceScale_(faceScale),
136  isBasisHCURL_(isBasisHCURL),
137  isBasisHDIV_(isBasisHDIV)
138  {
139  if(numEdges > 0)
140  edgeWorkView_ = ScalarViewType("edgeWorkView", dofCoords.extent(0), numEdgeInternalDofs, 1);
141  if(numFaces > 0)
142  faceWorkView_ = ScalarViewType("faceWorkView", dofCoords.extent(0), facesInternalDofCoords.extent(1), 2);
143  }
144 
145  KOKKOS_INLINE_FUNCTION
146  void operator()(const ordinal_type cell) const {
147  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
148 
149 
150  if(numEdges_ > 0) {
151  //compute coordinates associated to edge DoFs
152  ordinal_type eOrt[12];
153  value_type ortJac;
154  ScalarViewType ortJacobianEdge(&ortJac, 1, 1);
155  orts_(cell).getEdgeOrientation(eOrt, numEdges_);
156  auto edgeInternalDofCoordsOriented = Kokkos::subview(edgeWorkView_,cell, Kokkos::ALL(), Kokkos::ALL());
157  //map edge DoFs coords into parent element coords
158  for (ordinal_type iedge=0; iedge < numEdges_; ++iedge) {
159  Impl::OrientationTools::mapToModifiedReference(edgeInternalDofCoordsOriented,edgeInternalDofCoords_,edgeTopoKey_, eOrt[iedge]);
160 
161  for (ordinal_type j=0;j<numEdgeInternalDofs_;++j) {
162  const auto u = edgeInternalDofCoordsOriented(j, 0);
163  auto idof = tagToOrdinal_(1, iedge, j);
164  for (ordinal_type d=0;d<cellDim_;++d)
165  dofCoords_(cell,idof,d) = edgeParam_(iedge, d, 0) + edgeParam_(iedge, d, 1)*u ;
166  }
167  }
168 
169  //compute coefficients associated to edge DoFs
170  if(isBasisHCURL_) {
171  for (ordinal_type iedge=0; iedge < numEdges_; ++iedge) {
172  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobianEdge, edgeTopoKey_, eOrt[iedge]);
173  for(ordinal_type j=0; j<numEdgeInternalDofs_; ++j) {
174  auto idof = tagToOrdinal_(1, iedge, j);
175  auto jdof = edgeInternalDofOrdinals_(j);
176  for(ordinal_type d=0; d <cellDim_; ++d) {
177  dofCoeffs_(cell,idof,d) = 0;
178  for(ordinal_type k=0; k <1; ++k)
179  for(ordinal_type l=0; l <1; ++l)
180  dofCoeffs_(cell,idof,d) += refEdgesTan_(iedge,d)*ortJacobianEdge(0,0)*edgeDofCoeffs_(jdof)*edgeScale_;
181  }
182  }
183  }
184  } else if(isBasisHDIV_) {
185  for (ordinal_type iedge=0; iedge < numEdges_; ++iedge) {
186  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobianEdge, edgeTopoKey_, eOrt[iedge]);
187  auto ortJacobianDet = ortJacobianEdge(0,0);
188  for(ordinal_type j=0; j<numEdgeInternalDofs_; ++j) {
189  auto idof = tagToOrdinal_(1, iedge, j);
190  auto jdof = edgeInternalDofOrdinals_(j);
191  for(ordinal_type d=0; d <cellDim_; ++d)
192  dofCoeffs_(cell,idof,d) = refEdgesNormal_(iedge, d)*ortJacobianDet*edgeDofCoeffs_(jdof)*edgeScale_;
193  }
194  }
195  } else {
196  for (ordinal_type iedge=0; iedge < numEdges_; ++iedge) {
197  for(ordinal_type j=0; j<numEdgeInternalDofs_; ++j) {
198  auto idof = tagToOrdinal_(1, iedge, j);
199  auto jdof = edgeInternalDofOrdinals_(j);
200  dofCoeffs_(cell,idof,0) = edgeDofCoeffs_(jdof);
201  }
202  }
203 
204  }
205  }
206 
207  if(numFaces_ > 0) {
208  //compute coordinates associated to face DoFs
209  ordinal_type fOrt[12];
210  value_type ortJac[4];
211  ScalarViewType ortJacobianFace(ortJac, 2, 2);
212  orts_(cell).getFaceOrientation(fOrt, numFaces_);
213  //map face dofs coords into parent element coords
214  for (ordinal_type iface=0; iface < numFaces_; ++iface) {
215  ordinal_type ort = fOrt[iface];
216  ordinal_type numInternalDofs = numFacesInternalDofs_(iface);
217  auto dofRange = range_type(0, numInternalDofs);
218  auto faceInternalDofCoords = Kokkos::subview(facesInternalDofCoords_, iface, dofRange, Kokkos::ALL());
219  auto faceInternalDofCoordsOriented = Kokkos::subview(faceWorkView_,cell, dofRange, Kokkos::ALL());
220  Impl::OrientationTools::mapToModifiedReference(faceInternalDofCoordsOriented,faceInternalDofCoords,faceTopoKey_(iface),ort);
221 
222  for (ordinal_type j=0;j<numInternalDofs;++j) {
223  const auto u = faceInternalDofCoordsOriented(j, 0);
224  const auto v = faceInternalDofCoordsOriented(j, 1);
225  auto idof = tagToOrdinal_(2, iface, j);
226  for (ordinal_type d=0;d<cellDim_;++d)
227  dofCoords_(cell,idof,d) = faceParam_(iface, d, 0) + faceParam_(iface, d, 1)*u + faceParam_(iface, d, 2)*v;
228  }
229  }
230  //compute coefficients associated to face DoFs
231  if(isBasisHCURL_) {
232  for (ordinal_type iface=0; iface < numFaces_; ++iface) {
233  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobianFace, faceTopoKey_(iface), fOrt[iface]);
234  for(ordinal_type j=0; j<numFacesInternalDofs_(iface); ++j) {
235  auto idof = tagToOrdinal_(2, iface, j);
236  auto jdof = facesInternalDofOrdinals_(iface, j);
237  for(ordinal_type d=0; d <cellDim_; ++d) {
238  dofCoeffs_(cell,idof,d) = 0;
239  for(ordinal_type k=0; k <2; ++k)
240  for(ordinal_type l=0; l <2; ++l)
241  dofCoeffs_(cell,idof,d) += refFaceTangents_(iface, d,l)*ortJacobianFace(l,k)*faceDofCoeffs_(iface,jdof,k);
242  }
243  }
244  }
245  } else if(isBasisHDIV_) {
246  for (ordinal_type iface=0; iface < numFaces_; ++iface) {
247  Impl::OrientationTools::getJacobianOfOrientationMap(ortJacobianFace, faceTopoKey_(iface), fOrt[iface]);
248  auto ortJacobianDet = ortJacobianFace(0,0)*ortJacobianFace(1,1)-ortJacobianFace(1,0)*ortJacobianFace(0,1);
249  for(ordinal_type j=0; j<numFacesInternalDofs_(iface); ++j) {
250  auto idof = tagToOrdinal_(2, iface, j);
251  auto jdof = facesInternalDofOrdinals_(iface, j);
252  for(ordinal_type d=0; d <cellDim_; ++d)
253  dofCoeffs_(cell,idof,d) = refFacesNormal_(iface,d)*ortJacobianDet*faceDofCoeffs_(iface,jdof)*faceScale_;
254  }
255  }
256  } else {
257 
258  for (ordinal_type iface=0; iface < numFaces_; ++iface) {
259  for(ordinal_type j=0; j<numFacesInternalDofs_(iface); ++j) {
260  auto idof = tagToOrdinal_(2, iface, j);
261  auto jdof = facesInternalDofOrdinals_(iface, j);
262  dofCoeffs_(cell,idof,0) = faceDofCoeffs_(iface,jdof);
263  }
264  }
265  }
266  }
267  }
268 };
269 
270 template<typename SpT>
271 template<typename BasisType,
272 class ...coordsProperties, class ...coeffsProperties,
273 typename ortValueType, class ...ortProperties>
274 void
276  Kokkos::DynRankView<typename BasisType::scalarType, coordsProperties...> dofCoords,
277  Kokkos::DynRankView<typename BasisType::scalarType, coeffsProperties...> dofCoeffs,
278  const BasisType* basis,
279  EPointType pointType,
280  const Kokkos::DynRankView<ortValueType, ortProperties...> orts) {
281 
282  typedef typename BasisType::scalarType scalarType;
283  typedef Kokkos::DynRankView<scalarType, SpT> ScalarViewType;
284  typedef Kokkos::DynRankView<ordinal_type, SpT> intViewType;
285  typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
286 
287  const auto topo = basis->getBaseCellTopology();
288  const std::string name(basis->getName());
289  const ordinal_type degree(basis->getDegree());
290 
291  bool isBasisHCURL = name.find("HCURL") != std::string::npos;
292  bool isBasisHDIV = name.find("HDIV") != std::string::npos;
293  bool isBasisTriOrTet = (name.find("TRI") != std::string::npos) || (name.find("TET") != std::string::npos);
294  bool isBasisHEXI1 = (name.find("HEX_I1") != std::string::npos);
295  bool isBasisTETI1 = (name.find("TET_I1") != std::string::npos);
296  bool isBasisI1 = (name.find("I1") != std::string::npos);
297 
298  ordinal_type numEdges = (basis->getDofCount(1, 0) > 0) ? topo.getEdgeCount() : 0;
299  ordinal_type numFaces = (basis->getDofCount(2, 0) > 0) ? topo.getFaceCount() : 0;
300 
301  Teuchos::RCP<Basis<SpT,scalarType,scalarType> > faceBases[6];
302  Teuchos::RCP<Basis<SpT,scalarType,scalarType> > edgeBasis;
303 
304 
305  if ((name == "Intrepid2_HGRAD_QUAD_Cn_FEM") || (name == "Intrepid2_HGRAD_QUAD_C1_FEM")) {
306  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree, pointType) );
307  } else if ((name == "Intrepid2_HGRAD_HEX_Cn_FEM") || (name == "Intrepid2_HGRAD_HEX_C1_FEM")){
308  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree, pointType) );
309  faceBases[0] = Teuchos::rcp(new Basis_HGRAD_QUAD_Cn_FEM<SpT,scalarType,scalarType>(degree, pointType) );
310  for(ordinal_type i=1; i<6; ++i) faceBases[i]=faceBases[0];
311  } else if ((name == "Intrepid2_HGRAD_TRI_Cn_FEM") || (name == "Intrepid2_HGRAD_TRI_C1_FEM")) {
312  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree, pointType) );
313  } else if ((name == "Intrepid2_HGRAD_TET_Cn_FEM") || (name == "Intrepid2_HGRAD_TET_C1_FEM")) {
314  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree, pointType) );
315  faceBases[0] = Teuchos::rcp(new Basis_HGRAD_TRI_Cn_FEM<SpT,scalarType,scalarType>(degree, pointType) );
316  for(ordinal_type i=1; i<4; ++i) faceBases[i]=faceBases[0];
317  } else if ((name == "Intrepid2_HCURL_QUAD_In_FEM") || (name == "Intrepid2_HCURL_QUAD_I1_FEM")) {
318  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree-1, POINTTYPE_GAUSS) );
319  } else if ((name == "Intrepid2_HCURL_HEX_In_FEM") || (name == "Intrepid2_HCURL_HEX_I1_FEM")) {
320  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree-1, POINTTYPE_GAUSS) );
321  faceBases[0] = Teuchos::rcp(new Basis_HCURL_QUAD_In_FEM<SpT,scalarType,scalarType>(degree, pointType) );
322  for(ordinal_type i=1; i<6; ++i) faceBases[i]=faceBases[0];
323  } else if ((name == "Intrepid2_HCURL_TRI_In_FEM") || (name == "Intrepid2_HCURL_TRI_I1_FEM")) {
324  edgeBasis = Teuchos::rcp(new Basis_HVOL_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree-1, pointType) );
325  } else if ((name == "Intrepid2_HCURL_TET_In_FEM") || (name == "Intrepid2_HCURL_TET_I1_FEM")) {
326  edgeBasis = Teuchos::rcp(new Basis_HVOL_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree-1, pointType) );
327  faceBases[0] = Teuchos::rcp(new Basis_HCURL_TRI_In_FEM<SpT,scalarType,scalarType>(degree, pointType) );
328  for(ordinal_type i=1; i<4; ++i) faceBases[i]=faceBases[0];
329  } else if ((name == "Intrepid2_HDIV_QUAD_In_FEM") || (name == "Intrepid2_HDIV_QUAD_I1_FEM")) {
330  edgeBasis = Teuchos::rcp(new Basis_HGRAD_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree-1, POINTTYPE_GAUSS) );
331  } else if ((name == "Intrepid2_HDIV_HEX_In_FEM") || (name == "Intrepid2_HDIV_HEX_I1_FEM")) {
332  edgeBasis = Teuchos::null;
333  faceBases[0] = Teuchos::rcp(new Basis_HGRAD_QUAD_Cn_FEM<SpT,scalarType,scalarType>(degree-1, POINTTYPE_GAUSS) );
334  for(ordinal_type i=1; i<6; ++i) faceBases[i]=faceBases[0];
335  } else if ((name == "Intrepid2_HDIV_TRI_In_FEM") || (name == "Intrepid2_HDIV_TRI_I1_FEM")) {
336  edgeBasis = Teuchos::rcp(new Basis_HVOL_LINE_Cn_FEM<SpT,scalarType,scalarType>(degree-1, pointType) );
337  } else if ((name == "Intrepid2_HDIV_TET_In_FEM") || (name == "Intrepid2_HDIV_TET_I1_FEM")) {
338  edgeBasis = Teuchos::null;
339  faceBases[0] = Teuchos::rcp(new Basis_HVOL_TRI_Cn_FEM<SpT,scalarType,scalarType>(degree-1, pointType) );
340  for(ordinal_type i=1; i<4; ++i) faceBases[i]=faceBases[0];
341  } else { //HVOL element does not have any face or edge DOFs
342  //Throw error when basis is not HVOL.
343  INTREPID2_TEST_FOR_ABORT(name.find("HVOL") == std::string::npos,
344  ">>> ERROR (Intrepid2::Experimental::LagrangianInterpolation<SpT>::getDofCoordsAndCoeffs): " \
345  "method not implemented for this basis function");
346  }
347 
348  auto tagToOrdinal = Kokkos::create_mirror_view_and_copy(typename SpT::memory_space(), basis->getAllDofOrdinal());
349 
350  const ordinal_type dim = topo.getDimension();
351 
352  const ordinal_type numCells = dofCoeffs.extent(0);
353 
354  ScalarViewType refDofCoords("refDofCoords", dofCoords.extent(1), dofCoords.extent(2)), refDofCoeffs;
355  basis->getDofCoords(refDofCoords);
356  RealSpaceTools<SpT>::clone(dofCoords,refDofCoords);
357 
358  if(dofCoeffs.rank() == 3) //vector basis
359  refDofCoeffs = ScalarViewType("refDofCoeffs", dofCoeffs.extent(1), dofCoeffs.extent(2));
360  else //scalar basis
361  refDofCoeffs = ScalarViewType("refDofCoeffs",dofCoeffs.extent(1));
362  basis->getDofCoeffs(refDofCoeffs);
363  RealSpaceTools<SpT>::clone(dofCoeffs,refDofCoeffs);
364 
365  //*** Pre-compute needed quantities related to edge DoFs that do not depend on the cell ***
366 
367  ordinal_type edgeTopoKey = Teuchos::nonnull(edgeBasis) ? edgeBasis->getBaseCellTopology().getBaseKey() : 0;
368  intViewType eOrt("eOrt", numEdges);
369  ScalarViewType refEdgesTan("refEdgesTan", numEdges, dim);
370  ScalarViewType refEdgesNormal("refEdgesNormal", numEdges, dim);
371  ScalarViewType edgeParam;
372  ordinal_type edgeBasisCardinality = Teuchos::nonnull(edgeBasis) ? edgeBasis->getCardinality() : ordinal_type(0);
373  ordinal_type numEdgeInternalDofs = Teuchos::nonnull(edgeBasis) ? edgeBasis->getDofCount(1,0) : ordinal_type(0);
374  ScalarViewType edgeDofCoords("edgeDofCoords", edgeBasisCardinality, 1);
375  ScalarViewType edgeDofCoeffs("edgeDofCoeffs", edgeBasisCardinality);
376  ScalarViewType edgeInternalDofCoords("edgeInternalDofCoords", numEdgeInternalDofs, 1);
377  intViewType edgeInternalDofOrdinals("edgeInternalDofOrdinals", numEdgeInternalDofs);
378  //depending on how the reference basis is defined, the edges are scaled differently
379  auto edgeScale = (isBasisTriOrTet||isBasisI1) ? 2.0 :1.0;
380 
381  if(Teuchos::nonnull(edgeBasis)) {
382  edgeBasis->getDofCoords(edgeDofCoords);
383  edgeBasis->getDofCoeffs(edgeDofCoeffs);
384  }
385 
386  for (ordinal_type iedge=0; iedge < numEdges; ++iedge) {
387  if(isBasisHCURL) {
388  auto edgeTan = Kokkos::subview(refEdgesTan, iedge, Kokkos::ALL);
389  auto edgeTanHost = Kokkos::create_mirror_view(edgeTan);
390  CellTools<SpT>::getReferenceEdgeTangent(edgeTanHost, iedge, topo);
391  Kokkos::deep_copy(edgeTan,edgeTanHost);
392  }
393  else if(isBasisHDIV) {
394  auto edgeNormal = Kokkos::subview(refEdgesNormal, iedge, Kokkos::ALL);
395  auto edgeNormalHost = Kokkos::create_mirror_view(edgeNormal);
396  CellTools<SpT>::getReferenceSideNormal(edgeNormalHost, iedge, topo);
397  Kokkos::deep_copy(edgeNormal,edgeNormalHost);
398  }
399  }
400 
401  //compute DofCoords Oriented
402  for(ordinal_type i=0; i<numEdgeInternalDofs; ++i) {
403  edgeInternalDofOrdinals(i) = edgeBasis->getDofOrdinal(1, 0, i);
404  edgeInternalDofCoords(i,0) = edgeDofCoords(edgeInternalDofOrdinals(i), 0);
405  CellTools<SpT>::getSubcellParametrization(edgeParam, 1, topo);
406  }
407 
408 
409  //*** Pre-compute needed quantities related to face DoFs that do not depend on the cell ***
410 
411  intViewType faceTopoKey("faceTopoKey",numFaces);
412  intViewType fOrt("fOrt", numFaces);
413  ScalarViewType refFaceTangents("refFaceTangents", numFaces, dim, 2);
414  ScalarViewType refFacesNormal("refFacesNormal", numFaces, dim);
415  ScalarViewType faceParam;
416  intViewType numFacesInternalDofs("numFacesInternalDofs", numFaces);
417  ScalarViewType facesInternalDofCoords;
418  intViewType facesInternalDofOrdinals;
419  ScalarViewType faceDofCoeffs;
420  //depending on how the reference basis is defined, the faces are scaled differently
421  auto faceScale = (isBasisHEXI1) ? 4.0 :
422  (isBasisTETI1) ? 0.5 : 1.0;
423 
424 
425  ordinal_type maxNumFacesInternalDofs=0;
426  ordinal_type faceBasisMaxCardinality=0;
427 
428  for (ordinal_type iface=0; iface < numFaces; ++iface) {
429  ordinal_type numInternalDofs = faceBases[iface]->getDofCount(2,0);
430  numFacesInternalDofs(iface) = numInternalDofs;
431  maxNumFacesInternalDofs = std::max(maxNumFacesInternalDofs,numInternalDofs);
432  ordinal_type faceBasisCardinality = faceBases[iface]->getCardinality();
433  faceBasisMaxCardinality = std::max(faceBasisMaxCardinality, faceBasisCardinality);
434  }
435 
436  facesInternalDofCoords = ScalarViewType("faceInternalDofCoords", numFaces, maxNumFacesInternalDofs, 2);
437  facesInternalDofOrdinals = intViewType("faceInternalDofCoords", numFaces, maxNumFacesInternalDofs);
438 
439  if(isBasisHCURL)
440  faceDofCoeffs = ScalarViewType("faceDofCoeffs", numFaces, faceBasisMaxCardinality,2);
441  else
442  faceDofCoeffs = ScalarViewType("faceDofCoeffs", numFaces, faceBasisMaxCardinality);
443 
444  for (ordinal_type iface=0; iface < numFaces; ++iface) {
445  auto faceBasis = faceBases[iface];
446  faceTopoKey(iface) = faceBasis->getBaseCellTopology().getBaseKey();
447  ordinal_type faceBasisCardinality = faceBasis->getCardinality();
448  ScalarViewType faceDofCoords("faceDofCoords", faceBasisCardinality, 2);
449  faceBasis->getDofCoords(faceDofCoords);
450  for(ordinal_type i=0; i<numFacesInternalDofs(iface); ++i) {
451  facesInternalDofOrdinals(iface, i) = faceBasis->getDofOrdinal(2, 0, i);
452  for(ordinal_type d=0; d <2; ++d)
453  facesInternalDofCoords(iface, i,d) = faceDofCoords(facesInternalDofOrdinals(iface, i),d);
454  }
455 
456  auto dofRange = range_type(0, faceBasis->getCardinality());
457  faceBasis->getDofCoeffs(Kokkos::subview(faceDofCoeffs, iface, dofRange, Kokkos::ALL()));
458 
459  if(isBasisHCURL) {
460  auto refFaceTanU = Kokkos::subview(refFaceTangents, iface, Kokkos::ALL, 0);
461  auto refFaceTanV = Kokkos::subview(refFaceTangents, iface, Kokkos::ALL,1);
462  auto refFaceTanUHost = Kokkos::create_mirror_view(refFaceTanU);
463  auto refFaceTanVHost = Kokkos::create_mirror_view(refFaceTanV);
464  CellTools<SpT>::getReferenceFaceTangents(refFaceTanUHost, refFaceTanVHost, iface, topo);
465  Kokkos::deep_copy(refFaceTanU, refFaceTanUHost);
466  Kokkos::deep_copy(refFaceTanV, refFaceTanVHost);
467  } else if(isBasisHDIV) {
468  auto faceNormal = Kokkos::subview(refFacesNormal,iface,Kokkos::ALL());
469  auto faceNormalHost = Kokkos::create_mirror_view(faceNormal);
470  CellTools<SpT>::getReferenceFaceNormal(faceNormalHost, iface, topo);
471  Kokkos::deep_copy(faceNormal, faceNormalHost);
472  }
473  }
474 
475  if(dim > 2)
476  CellTools<SpT>::getSubcellParametrization(faceParam, 2, topo);
477 
478 
479  //*** Loop over cells ***
480 
481  const Kokkos::RangePolicy<SpT> policy(0, numCells);
483  <ScalarViewType,
484  decltype(orts),
485  decltype(tagToOrdinal),
486  decltype(edgeParam),
487  intViewType> FunctorType;
488  Kokkos::parallel_for(policy,
489  FunctorType(dofCoords, dofCoeffs,
490  orts, tagToOrdinal, edgeParam, faceParam,
491  edgeInternalDofOrdinals, facesInternalDofOrdinals,
492  edgeInternalDofCoords, edgeDofCoeffs, refEdgesTan, refEdgesNormal,
493  facesInternalDofCoords, faceDofCoeffs, refFaceTangents, refFacesNormal,
494  dim, numEdges, numFaces,
495  edgeTopoKey, numEdgeInternalDofs,
496  faceTopoKey, numFacesInternalDofs,
497  edgeScale, faceScale,
498  isBasisHCURL, isBasisHDIV));
499 }
500 
501 
502 template<typename SpT>
503 template<typename basisCoeffsViewType,
504 typename funcViewType,
505 typename dofCoeffViewType>
506 void
507 LagrangianInterpolation<SpT>::getBasisCoeffs(basisCoeffsViewType basisCoeffs,
508  const funcViewType functionValsAtDofCoords,
509  const dofCoeffViewType dofCoeffs){
510  ArrayTools<SpT>::dotMultiplyDataData(basisCoeffs,functionValsAtDofCoords,dofCoeffs);
511 }
512 }
513 }
514 
515 #endif
516 
Implementation of the default HVOL-compatible Lagrange basis of arbitrary degree on Triangle cell...
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.
Implementation of the default H(curl)-compatible FEM basis on Quadrilateral cell. ...
static void getReferenceSideNormal(Kokkos::DynRankView< refSideNormalValueType, refSideNormalProperties...> refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
static void getDofCoordsAndCoeffs(Kokkos::DynRankView< typename BasisType::scalarType, coordsProperties...> dofCoords, Kokkos::DynRankView< typename BasisType::scalarType, coeffsProperties...> dofCoeffs, const BasisType *cellBasis, EPointType basisPointType, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations)
Computes the points and coefficients associated with the basis DOFs for the reference oriented elemen...
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 dotMultiplyDataData(Kokkos::DynRankView< outputDataValueType, outputDataProperties...> outputData, const Kokkos::DynRankView< inputDataLeftValueType, inputDataLeftProperties...> inputDataLeft, const Kokkos::DynRankView< inputDataRightValueType, inputDataRightProperties...> inputDataRight)
There are two use cases: (1) dot product of a rank-2, 3 or 4 container inputDataRight with dimensions...
Header file for the Intrepid2::FunctionSpaceTools class.
static void clone(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input)
Clone input array.
Implementation of the locally H(grad)-compatible FEM basis of variable order on the [-1...
Header file for Intrepid2::ArrayTools class providing utilities for array operations.
Implementation of the locally HVOL-compatible FEM basis of variable order on the [-1,1] reference line cell, using Lagrange polynomials.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Tr...
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
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 getBasisCoeffs(basisCoeffsViewType basisCoeffs, const funcViewType functionAtDofCoords, const dofCoeffViewType dofCoeffs)
Computes the basis weights of the function interpolation.
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...
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.
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...