Intrepid2
Intrepid2_OrientationToolsDefMatrixData.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 
15 #ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_MATRIX_DATA_HPP__
16 #define __INTREPID2_ORIENTATIONTOOLS_DEF_MATRIX_DATA_HPP__
17 
18 // disable clang warnings
19 #if defined (__clang__) && !defined (__INTEL_COMPILER)
20 #pragma clang system_header
21 #endif
22 
23 namespace Intrepid2 {
24 
25  template<typename DT>
26  template<typename BasisHostType>
27  typename OrientationTools<DT>::CoeffMatrixDataViewType
28  OrientationTools<DT>::createCoeffMatrixInternal(const BasisHostType* basis, const bool inverse) {
29  const std::string name(basis->getName());
30  CoeffMatrixDataViewType matData;
31 
32  const auto cellTopo = basis->getBaseCellTopology();
33  const ordinal_type numEdges = cellTopo.getSubcellCount(1);
34  const ordinal_type numFaces = cellTopo.getSubcellCount(2);
35  ordinal_type matDim = 0, matDim1 = 0, matDim2 = 0, numOrts = 0, numSubCells;
36  for(ordinal_type i=0; i<numEdges; ++i) {
37  matDim1 = std::max(matDim1, basis->getDofCount(1,i));
38  numOrts = std::max(numOrts,2);
39  }
40  for(ordinal_type i=0; i<numFaces; ++i) {
41  matDim2 = std::max(matDim2, basis->getDofCount(2,i));
42  numOrts = std::max(numOrts,2*ordinal_type(cellTopo.getSideCount(2,i)));
43  }
44  matDim = std::max(matDim1,matDim2);
45  numSubCells = (matDim1>0)*numEdges + (matDim2>0)*numFaces;
46 
47 
48  matData = CoeffMatrixDataViewType("Orientation::CoeffMatrix::"+name,
49  numSubCells,
50  numOrts,
51  matDim,
52  matDim);
53 
54  if(basis->getFunctionSpace() == FUNCTION_SPACE_HGRAD) {
55  init_HGRAD(matData, basis, inverse);
56  } else if (basis->getFunctionSpace() == FUNCTION_SPACE_HCURL) {
57  init_HCURL(matData, basis, inverse);
58  } else if (basis->getFunctionSpace() == FUNCTION_SPACE_HDIV) {
59  init_HDIV(matData, basis, inverse);
60  } else if (basis->getFunctionSpace() == FUNCTION_SPACE_HVOL) {
61  init_HVOL(matData, basis, inverse);
62  }
63  return matData;
64  }
65 
66  //
67  // HGRAD elements
68  //
69  template<typename DT>
70  template<typename BasisHostType>
71  void
73  init_HGRAD(typename OrientationTools<DT>::CoeffMatrixDataViewType matData,
74  BasisHostType const *cellBasis,
75  const bool inverse) {
76 
77  const auto cellTopo = cellBasis->getBaseCellTopology();
78  const ordinal_type numEdges = cellTopo.getSubcellCount(1);
79  const ordinal_type numFaces = cellTopo.getSubcellCount(2);
80  Intrepid2::BasisPtr<typename BasisHostType::DeviceType,
81  typename BasisHostType::OutputValueType,
82  typename BasisHostType::PointValueType>
83  basisPtr;
84  BasisHostType const *subcellBasis;
85 
86  { //edges
87  subcellBasis = cellBasis; // if (dim==1)
88  const ordinal_type numOrt = 2;
89  for (ordinal_type edgeId=0;edgeId<numEdges;++edgeId) {
90  if(cellBasis->getDofCount(1, edgeId) < 2) continue;
91  if(cellTopo.getDimension()!=1) {
92  basisPtr = cellBasis->getSubCellRefBasis(1,edgeId);
93  subcellBasis = basisPtr.get();
94  }
95 
96  for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
97  auto mat = Kokkos::subview(matData,
98  edgeId, edgeOrt,
99  Kokkos::ALL(), Kokkos::ALL());
101  (mat,
102  *subcellBasis, *cellBasis,
103  edgeId, edgeOrt, inverse);
104  }
105  }
106  }
107  { //faces
108  subcellBasis = cellBasis; // if(dim==2)
109  for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
110  // this works for triangles (numOrt=6) and quadrilaterals (numOrt=8)
111  const ordinal_type numOrt = 2*cellTopo.getSideCount(2,faceId);
112  if(cellBasis->getDofCount(2, faceId) < 1) continue;
113  if(cellTopo.getDimension()!=2) {
114  basisPtr = cellBasis->getSubCellRefBasis(2,faceId);
115  subcellBasis = basisPtr.get();
116  }
117  for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
118  auto mat = Kokkos::subview(matData,
119  numEdges+faceId, faceOrt,
120  Kokkos::ALL(), Kokkos::ALL());
122  (mat,
123  *subcellBasis, *cellBasis,
124  faceId, faceOrt, inverse);
125  }
126  }
127  }
128  }
129 
130  //
131  // HCURL elements
132  //
133  template<typename DT>
134  template<typename BasisHostType>
135  void
137  init_HCURL(typename OrientationTools<DT>::CoeffMatrixDataViewType matData,
138  BasisHostType const *cellBasis, const bool inverse) {
139  const auto cellTopo = cellBasis->getBaseCellTopology();
140  const ordinal_type numEdges = cellTopo.getSubcellCount(1);
141  const ordinal_type numFaces = cellTopo.getSubcellCount(2);
142  Intrepid2::BasisPtr<typename BasisHostType::DeviceType,
143  typename BasisHostType::OutputValueType,
144  typename BasisHostType::PointValueType>
145  basisPtr;
146  BasisHostType const* subcellBasis;
147 
148  { // edges
149  subcellBasis = cellBasis; // if (dim==1)
150  const ordinal_type numOrt = 2;
151  for (ordinal_type edgeId=0;edgeId<numEdges;++edgeId) {
152  if(cellBasis->getDofCount(1, edgeId) < 1) continue;
153  if(cellTopo.getDimension()!=1) {
154  basisPtr = cellBasis->getSubCellRefBasis(1,edgeId);
155  subcellBasis = basisPtr.get();
156  }
157  for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
158  auto mat = Kokkos::subview(matData,
159  edgeId, edgeOrt,
160  Kokkos::ALL(), Kokkos::ALL());
162  *subcellBasis, *cellBasis,
163  edgeId, edgeOrt, inverse);
164  }
165  }
166  }
167  { //faces
168  subcellBasis = cellBasis; // if (dim==2)
169  for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
170  // this works for triangles (numOrt=6) and quadratures (numOrt=8)
171  const ordinal_type numOrt = 2*cellTopo.getSideCount(2,faceId);
172  if(cellBasis->getDofCount(2, faceId) < 1) continue;
173  if(cellTopo.getDimension()!=2) {
174  basisPtr = cellBasis->getSubCellRefBasis(2,faceId);
175  subcellBasis = basisPtr.get();
176  }
177  for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
178  auto mat = Kokkos::subview(matData,
179  numEdges+faceId, faceOrt,
180  Kokkos::ALL(), Kokkos::ALL());
182  (mat,
183  *subcellBasis, *cellBasis,
184  faceId, faceOrt, inverse);
185  }
186  }
187  }
188  }
189 
190  //
191  // HDIV elements
192  //
193  template<typename DT>
194  template<typename BasisHostType>
195  void
197  init_HDIV(typename OrientationTools<DT>::CoeffMatrixDataViewType matData,
198  BasisHostType const *cellBasis, const bool inverse) {
199  const auto cellTopo = cellBasis->getBaseCellTopology();
200  const ordinal_type numSides = cellTopo.getSideCount();
201  const ordinal_type sideDim = cellTopo.getDimension()-1;
202  Intrepid2::BasisPtr<typename BasisHostType::DeviceType,
203  typename BasisHostType::OutputValueType,
204  typename BasisHostType::PointValueType>
205  subcellBasisPtr;
206 
207  {
208  for (ordinal_type sideId=0;sideId<numSides;++sideId) {
209  if(cellBasis->getDofCount(sideDim, sideId) < 1) continue;
210  const ordinal_type numOrt = (sideDim == 1) ? 2 : 2*cellTopo.getSideCount(sideDim,sideId);
211  subcellBasisPtr = cellBasis->getSubCellRefBasis(sideDim,sideId);
212  for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
213  auto mat = Kokkos::subview(matData,
214  sideId, faceOrt,
215  Kokkos::ALL(), Kokkos::ALL());
217  *subcellBasisPtr, *cellBasis,
218  sideId, faceOrt, inverse);
219  }
220  }
221  }
222  }
223 
224  //
225  // HVOL elements (used for 2D and 1D side cells)
226  //
227  template<typename DT>
228  template<typename BasisHostType>
229  void
231  init_HVOL(typename OrientationTools<DT>::CoeffMatrixDataViewType matData,
232  BasisHostType const *cellBasis, const bool inverse) {
233 
234  const auto cellTopo = cellBasis->getBaseCellTopology();
235  const ordinal_type numEdges = (cellTopo.getDimension()==1);
236  const ordinal_type numFaces = (cellTopo.getDimension()==2);
237 
238  {
239  const ordinal_type numOrt = 2;
240  for (ordinal_type edgeId=0;edgeId<numEdges;++edgeId) {
241  if(cellBasis->getDofCount(1, edgeId) < 1) continue;
242  for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
243  auto mat = Kokkos::subview(matData,
244  edgeId, edgeOrt,
245  Kokkos::ALL(), Kokkos::ALL());
247  (mat, *cellBasis, edgeOrt, inverse);
248  }
249  }
250  }
251  {
252  for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
253  // this works for triangles (numOrt=6) and quadratures (numOrt=8)
254  const ordinal_type numOrt = 2*cellTopo.getSideCount(2,faceId);
255  if(cellBasis->getDofCount(2, faceId) < 1) continue;
256  for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
257  auto mat = Kokkos::subview(matData,
258  numEdges+faceId, faceOrt,
259  Kokkos::ALL(), Kokkos::ALL());
261  (mat, *cellBasis, faceOrt, inverse);
262  }
263  }
264  }
265  }
266 
267  template<typename DT>
268  template<typename BasisType>
269  typename OrientationTools<DT>::CoeffMatrixDataViewType
270  OrientationTools<DT>::createCoeffMatrix(const BasisType* basis) {
271  Kokkos::push_finalize_hook( [=] {
272  ortCoeffData=std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<DT>::CoeffMatrixDataViewType>();
273  });
274 
275  const std::pair<std::string,ordinal_type> key(basis->getName(), basis->getDegree());
276  const auto found = ortCoeffData.find(key);
277 
278  CoeffMatrixDataViewType matData;
279  if (found == ortCoeffData.end()) {
280  {
281  auto basis_host = basis->getHostBasis();
282  matData = createCoeffMatrixInternal(basis_host.getRawPtr());
283  ortCoeffData.insert(std::make_pair(key, matData));
284  }
285  } else {
286  matData = found->second;
287  }
288 
289  return matData;
290  }
291 
292  template<typename DT>
293  template<typename BasisType>
294  typename OrientationTools<DT>::CoeffMatrixDataViewType
295  OrientationTools<DT>::createInvCoeffMatrix(const BasisType* basis) {
296  Kokkos::push_finalize_hook( [=] {
297  ortInvCoeffData=std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<DT>::CoeffMatrixDataViewType>();
298  });
299 
300  const std::pair<std::string,ordinal_type> key(basis->getName(), basis->getDegree());
301  const auto found = ortInvCoeffData.find(key);
302 
303  CoeffMatrixDataViewType matData;
304  if (found == ortInvCoeffData.end()) {
305  {
306  auto basis_host = basis->getHostBasis();
307  matData = createCoeffMatrixInternal(basis_host.getRawPtr(),true);
308  ortInvCoeffData.insert(std::make_pair(key, matData));
309  }
310  } else {
311  matData = found->second;
312  }
313 
314  return matData;
315  }
316 
317  template<typename DT>
319  ortCoeffData.clear();
320  ortInvCoeffData.clear();
321  }
322 }
323 
324 #endif
static void getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HGRAD basis for a given subcell and its reference basis.
static void getCoeffMatrix_HCURL(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HCURL basis for a given subcell and its reference basis.
static void clearCoeffMatrix()
Clear coefficient matrix.
static void init_HDIV(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HDIV basis.
static void getCoeffMatrix_HVOL(OutputViewType &output, const cellBasisHostType &cellBasis, const ordinal_type cellOrt, const bool inverse=false)
Compute orientation matrix for HVOL basis for a given (2D or 1D) cell and its reference basis...
static void init_HGRAD(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HGRAD basis.
static void getCoeffMatrix_HDIV(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
Compute orientation matrix for HDIV basis for a given subcell and its reference basis.
static CoeffMatrixDataViewType createCoeffMatrix(const BasisType *basis)
Create coefficient matrix.
static void init_HCURL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HCURL basis.
static void init_HVOL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
Compute orientation matrix for HVOL basis.
static CoeffMatrixDataViewType createInvCoeffMatrix(const BasisType *basis)
Create inverse of coefficient matrix.