Intrepid2
Intrepid2_OrientationToolsDefMatrixData.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
43
48#ifndef __INTREPID2_ORIENTATIONTOOLS_DEF_MATRIX_DATA_HPP__
49#define __INTREPID2_ORIENTATIONTOOLS_DEF_MATRIX_DATA_HPP__
50
51// disable clang warnings
52#if defined (__clang__) && !defined (__INTEL_COMPILER)
53#pragma clang system_header
54#endif
55
56namespace Intrepid2 {
57
58 template<typename DT>
59 template<typename BasisHostType>
61 OrientationTools<DT>::createCoeffMatrixInternal(const BasisHostType* basis) {
62 const std::string name(basis->getName());
63 CoeffMatrixDataViewType matData;
64
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);
72 }
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)));
76 }
77 matDim = std::max(matDim1,matDim2);
78 numSubCells = (matDim1>0)*numEdges + (matDim2>0)*numFaces;
79
80
81 matData = CoeffMatrixDataViewType("Orientation::CoeffMatrix::"+name,
82 numSubCells,
83 numOrts,
84 matDim,
85 matDim);
86
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);
95 }
96 return matData;
97 }
98
99 //
100 // HGRAD elements
101 //
102 template<typename DT>
103 template<typename BasisHostType>
104 void
107 BasisHostType const *cellBasis) {
108
109 const auto cellTopo = cellBasis->getBaseCellTopology();
110 const ordinal_type numEdges = cellTopo.getSubcellCount(1);
111 const ordinal_type numFaces = cellTopo.getSubcellCount(2);
112 Intrepid2::BasisPtr<typename BasisHostType::DeviceType,
113 typename BasisHostType::OutputValueType,
114 typename BasisHostType::PointValueType>
115 basisPtr;
116 BasisHostType const *subcellBasis;
117
118 { //edges
119 subcellBasis = cellBasis; // if (dim==1)
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();
126 }
127
128 for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
129 auto mat = Kokkos::subview(matData,
130 edgeId, edgeOrt,
131 Kokkos::ALL(), Kokkos::ALL());
133 (mat,
134 *subcellBasis, *cellBasis,
135 edgeId, edgeOrt);
136 }
137 }
138 }
139 { //faces
140 subcellBasis = cellBasis; // if(dim==2)
141 for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
142 // this works for triangles (numOrt=6) and quadrilaterals (numOrt=8)
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();
148 }
149 for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
150 auto mat = Kokkos::subview(matData,
151 numEdges+faceId, faceOrt,
152 Kokkos::ALL(), Kokkos::ALL());
154 (mat,
155 *subcellBasis, *cellBasis,
156 faceId, faceOrt);
157 }
158 }
159 }
160 }
161
162 //
163 // HCURL elements
164 //
165 template<typename DT>
166 template<typename BasisHostType>
167 void
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);
174 Intrepid2::BasisPtr<typename BasisHostType::DeviceType,
175 typename BasisHostType::OutputValueType,
176 typename BasisHostType::PointValueType>
177 basisPtr;
178 BasisHostType const* subcellBasis;
179
180 { // edges
181 subcellBasis = cellBasis; // if (dim==1)
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();
188 }
189 for (ordinal_type edgeOrt=0;edgeOrt<numOrt;++edgeOrt) {
190 auto mat = Kokkos::subview(matData,
191 edgeId, edgeOrt,
192 Kokkos::ALL(), Kokkos::ALL());
194 *subcellBasis, *cellBasis,
195 edgeId, edgeOrt);
196 }
197 }
198 }
199 { //faces
200 subcellBasis = cellBasis; // if (dim==2)
201 for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
202 // this works for triangles (numOrt=6) and quadratures (numOrt=8)
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();
208 }
209 for (ordinal_type faceOrt=0;faceOrt<numOrt;++faceOrt) {
210 auto mat = Kokkos::subview(matData,
211 numEdges+faceId, faceOrt,
212 Kokkos::ALL(), Kokkos::ALL());
214 (mat,
215 *subcellBasis, *cellBasis,
216 faceId, faceOrt);
217 }
218 }
219 }
220 }
221
222 //
223 // HDIV elements
224 //
225 template<typename DT>
226 template<typename BasisHostType>
227 void
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;
234 Intrepid2::BasisPtr<typename BasisHostType::DeviceType,
235 typename BasisHostType::OutputValueType,
236 typename BasisHostType::PointValueType>
237 subcellBasisPtr;
238
239 {
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,
246 sideId, faceOrt,
247 Kokkos::ALL(), Kokkos::ALL());
249 *subcellBasisPtr, *cellBasis,
250 sideId, faceOrt);
251 }
252 }
253 }
254 }
255
256 //
257 // HVOL elements (used for 2D and 1D side cells)
258 //
259 template<typename DT>
260 template<typename BasisHostType>
261 void
264 BasisHostType const *cellBasis) {
265
266 const auto cellTopo = cellBasis->getBaseCellTopology();
267 const ordinal_type numEdges = (cellTopo.getDimension()==1);
268 const ordinal_type numFaces = (cellTopo.getDimension()==2);
269
270 {
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,
276 edgeId, edgeOrt,
277 Kokkos::ALL(), Kokkos::ALL());
279 (mat, *cellBasis, edgeOrt);
280 }
281 }
282 }
283 {
284 for (ordinal_type faceId=0;faceId<numFaces;++faceId) {
285 // this works for triangles (numOrt=6) and quadratures (numOrt=8)
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);
294 }
295 }
296 }
297 }
298
299 template<typename DT>
300 template<typename BasisType>
303 Kokkos::push_finalize_hook( [=] {
304 ortCoeffData=std::map<std::pair<std::string,ordinal_type>, typename OrientationTools<DT>::CoeffMatrixDataViewType>();
305 });
306
307 const std::pair<std::string,ordinal_type> key(basis->getName(), basis->getDegree());
308 const auto found = ortCoeffData.find(key);
309
311 if (found == ortCoeffData.end()) {
312 {
313 auto basis_host = basis->getHostBasis();
314 matData = createCoeffMatrixInternal(basis_host.getRawPtr());
315 ortCoeffData.insert(std::make_pair(key, matData));
316 }
317 } else {
318 matData = found->second;
319 }
320
321 return matData;
322 }
323
324 template<typename DT>
328}
329
330#endif
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
static void getCoeffMatrix_HGRAD(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute orientation matrix for HGRAD basis for a given subcell and its reference basis.
static void getCoeffMatrix_HVOL(OutputViewType &output, const cellBasisHostType &cellBasis, const ordinal_type cellOrt)
Compute orientation matrix for HVOL basis for a given (2D or 1D) cell and its reference basis....
static void getCoeffMatrix_HCURL(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
Compute orientation matrix for HCURL basis for a given subcell and its reference basis.
static void getCoeffMatrix_HDIV(OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
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_HDIV(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HDIV basis.
static std::map< std::pair< std::string, ordinal_type >, CoeffMatrixDataViewType > ortCoeffData
key :: basis name, order, value :: matrix data view type
Kokkos::View< double ****, DeviceType > CoeffMatrixDataViewType
subcell ordinal, orientation, matrix m x n
static void init_HVOL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HVOL basis.
static void clearCoeffMatrix()
Clear coefficient matrix.
static void init_HCURL(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HCURL basis.
static void init_HGRAD(CoeffMatrixDataViewType matData, BasisHostType const *cellBasis)
Compute orientation matrix for HGRAD basis.