48#ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_MATRIX_DATA_HPP__
49#define __INTREPID2_ORIENTATIONTOOLS_DEF_MATRIX_DATA_HPP__
52#if defined (__clang__) && !defined (__INTEL_COMPILER)
53#pragma clang system_header
59 template<
typename BasisHostType>
61 OrientationTools<DT>::createCoeffMatrixInternal(
const BasisHostType* basis) {
62 const std::string name(basis->getName());
63 CoeffMatrixDataViewType matData;
65 const auto cellTopo = basis->getBaseCellTopology();
66 const ordinal_type numEdges = cellTopo.getSubcellCount(1);
67 const ordinal_type numFaces = cellTopo.getSubcellCount(2);
68 ordinal_type matDim = 0, matDim1 = 0, matDim2 = 0, numOrts = 0, numSubCells;
69 for(ordinal_type i=0; i<numEdges; ++i) {
70 matDim1 = std::max(matDim1, basis->getDofCount(1,i));
71 numOrts = std::max(numOrts,2);
73 for(ordinal_type i=0; i<numFaces; ++i) {
74 matDim2 = std::max(matDim2, basis->getDofCount(2,i));
75 numOrts = std::max(numOrts,2*ordinal_type(cellTopo.getSideCount(2,i)));
77 matDim = std::max(matDim1,matDim2);
78 numSubCells = (matDim1>0)*numEdges + (matDim2>0)*numFaces;
81 matData = CoeffMatrixDataViewType(
"Orientation::CoeffMatrix::"+name,
87 if(basis->getFunctionSpace() == FUNCTION_SPACE_HGRAD) {
88 init_HGRAD(matData, basis);
89 }
else if (basis->getFunctionSpace() == FUNCTION_SPACE_HCURL) {
90 init_HCURL(matData, basis);
91 }
else if (basis->getFunctionSpace() == FUNCTION_SPACE_HDIV) {
92 init_HDIV(matData, basis);
93 }
else if (basis->getFunctionSpace() == FUNCTION_SPACE_HVOL) {
94 init_HVOL(matData, basis);
102 template<
typename DT>
103 template<
typename BasisHostType>
107 BasisHostType
const *cellBasis) {
109 const auto cellTopo = cellBasis->getBaseCellTopology();
110 const ordinal_type numEdges = cellTopo.getSubcellCount(1);
111 const ordinal_type numFaces = cellTopo.getSubcellCount(2);
113 typename BasisHostType::OutputValueType,
114 typename BasisHostType::PointValueType>
116 BasisHostType
const *subcellBasis;
119 subcellBasis = cellBasis;
120 const ordinal_type numOrt = 2;
121 for (ordinal_type edgeId=0;edgeId<numEdges;++edgeId) {
122 if(cellBasis->getDofCount(1, edgeId) < 2)
continue;
123 if(cellTopo.getDimension()!=1) {
124 basisPtr = cellBasis->getSubCellRefBasis(1,edgeId);
125 subcellBasis = basisPtr.get();
128 for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
129 auto mat = Kokkos::subview(matData,
131 Kokkos::ALL(), Kokkos::ALL());
134 *subcellBasis, *cellBasis,
140 subcellBasis = cellBasis;
141 for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
143 const ordinal_type numOrt = 2*cellTopo.getSideCount(2,faceId);
144 if(cellBasis->getDofCount(2, faceId) < 1)
continue;
145 if(cellTopo.getDimension()!=2) {
146 basisPtr = cellBasis->getSubCellRefBasis(2,faceId);
147 subcellBasis = basisPtr.get();
149 for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
150 auto mat = Kokkos::subview(matData,
151 numEdges+faceId, faceOrt,
152 Kokkos::ALL(), Kokkos::ALL());
155 *subcellBasis, *cellBasis,
165 template<
typename DT>
166 template<
typename BasisHostType>
170 BasisHostType
const *cellBasis) {
171 const auto cellTopo = cellBasis->getBaseCellTopology();
172 const ordinal_type numEdges = cellTopo.getSubcellCount(1);
173 const ordinal_type numFaces = cellTopo.getSubcellCount(2);
175 typename BasisHostType::OutputValueType,
176 typename BasisHostType::PointValueType>
178 BasisHostType
const* subcellBasis;
181 subcellBasis = cellBasis;
182 const ordinal_type numOrt = 2;
183 for (ordinal_type edgeId=0;edgeId<numEdges;++edgeId) {
184 if(cellBasis->getDofCount(1, edgeId) < 1)
continue;
185 if(cellTopo.getDimension()!=1) {
186 basisPtr = cellBasis->getSubCellRefBasis(1,edgeId);
187 subcellBasis = basisPtr.get();
189 for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
190 auto mat = Kokkos::subview(matData,
192 Kokkos::ALL(), Kokkos::ALL());
194 *subcellBasis, *cellBasis,
200 subcellBasis = cellBasis;
201 for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
203 const ordinal_type numOrt = 2*cellTopo.getSideCount(2,faceId);
204 if(cellBasis->getDofCount(2, faceId) < 1)
continue;
205 if(cellTopo.getDimension()!=2) {
206 basisPtr = cellBasis->getSubCellRefBasis(2,faceId);
207 subcellBasis = basisPtr.get();
209 for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
210 auto mat = Kokkos::subview(matData,
211 numEdges+faceId, faceOrt,
212 Kokkos::ALL(), Kokkos::ALL());
215 *subcellBasis, *cellBasis,
225 template<
typename DT>
226 template<
typename BasisHostType>
230 BasisHostType
const *cellBasis) {
231 const auto cellTopo = cellBasis->getBaseCellTopology();
232 const ordinal_type numSides = cellTopo.getSideCount();
233 const ordinal_type sideDim = cellTopo.getDimension()-1;
235 typename BasisHostType::OutputValueType,
236 typename BasisHostType::PointValueType>
240 for (ordinal_type sideId=0;sideId<numSides;++sideId) {
241 if(cellBasis->getDofCount(sideDim, sideId) < 1)
continue;
242 const ordinal_type numOrt = (sideDim == 1) ? 2 : 2*cellTopo.getSideCount(sideDim,sideId);
243 subcellBasisPtr = cellBasis->getSubCellRefBasis(sideDim,sideId);
244 for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
245 auto mat = Kokkos::subview(matData,
247 Kokkos::ALL(), Kokkos::ALL());
249 *subcellBasisPtr, *cellBasis,
259 template<
typename DT>
260 template<
typename BasisHostType>
264 BasisHostType
const *cellBasis) {
266 const auto cellTopo = cellBasis->getBaseCellTopology();
267 const ordinal_type numEdges = (cellTopo.getDimension()==1);
268 const ordinal_type numFaces = (cellTopo.getDimension()==2);
271 const ordinal_type numOrt = 2;
272 for (ordinal_type edgeId=0;edgeId<numEdges;++edgeId) {
273 if(cellBasis->getDofCount(1, edgeId) < 1)
continue;
274 for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
275 auto mat = Kokkos::subview(matData,
277 Kokkos::ALL(), Kokkos::ALL());
279 (mat, *cellBasis, edgeOrt);
284 for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
286 const ordinal_type numOrt = 2*cellTopo.getSideCount(2,faceId);
287 if(cellBasis->getDofCount(2, faceId) < 1)
continue;
288 for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
289 auto mat = Kokkos::subview(matData,
290 numEdges+faceId, faceOrt,
291 Kokkos::ALL(), Kokkos::ALL());
293 (mat, *cellBasis, faceOrt);
299 template<
typename DT>
300 template<
typename BasisType>
303 Kokkos::push_finalize_hook( [=] {
307 const std::pair<std::string,ordinal_type> key(basis->getName(), basis->getDegree());
313 auto basis_host = basis->getHostBasis();
314 matData = createCoeffMatrixInternal(basis_host.getRawPtr());
318 matData = found->second;
324 template<
typename DT>
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.