Intrepid2
Intrepid2_ProjectionStructDef.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 
47 #ifndef __INTREPID2_PROJECTIONSTRUCTDEF_HPP__
48 #define __INTREPID2_PROJECTIONSTRUCTDEF_HPP__
49 
50 #include "Intrepid2_Utils.hpp"
51 #include "Shards_CellTopology.hpp"
52 #include "Shards_BasicTopologies.hpp"
53 
58 
60 
61 
62 namespace Intrepid2 {
63 
64 namespace Experimental {
65 template<typename SpT, typename ValueType>
66 template<typename BasisPtrType>
68  const ordinal_type targetCubDegree) {
69  const auto cellTopo = cellBasis->getBaseCellTopology();
70  std::string name = cellBasis->getName();
71  ordinal_type dim = cellTopo.getDimension();
72  numBasisEvalPoints = 0;
73  numBasisDerivEvalPoints = 0;
74  numTargetEvalPoints = 0;
75  numTargetDerivEvalPoints = 0;
76  const ordinal_type edgeDim = 1;
77  const ordinal_type faceDim = 2;
78 
79  ordinal_type basisCubDegree = cellBasis->getDegree();
80  ordinal_type edgeBasisCubDegree, faceBasisCubDegree;
81 
82  ordinal_type numVertices = (cellBasis->getDofCount(0, 0) > 0) ? cellTopo.getVertexCount() : 0;
83  ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount() : 0;
84  ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;
85 
86  ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);
87 
88  INTREPID2_TEST_FOR_ABORT( (numVertices > maxSubCellsCount) || (numFaces > maxSubCellsCount) || (numEdges > maxSubCellsCount),
89  ">>> ERROR (Intrepid2::ProjectionStruct:createHDivProjectionStruct, Projections do not support a cell topology with this cub cells count");
90 
91 
92  numBasisEvalPoints += numVertices;
93  numTargetEvalPoints += numVertices;
94  view_type coord("vertex_coord", dim);
95 
96  basisPointsRange = range_tag("basisPointsRange", 4,maxSubCellsCount);
97  targetPointsRange = range_tag("targetPointsRange", 4,maxSubCellsCount);
98  basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,maxSubCellsCount);
99  targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
100  subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);
101 
102  maxNumBasisEvalPoints = numVertices; maxNumTargetEvalPoints = numVertices;
103  for(ordinal_type iv=0; iv<numVertices; ++iv) {
104  basisPointsRange(0,iv) = range_type(iv, iv+1);
105  basisCubPoints[0][iv] = view_type("basisCubPoints",1,dim);
106  targetPointsRange(0,iv) = range_type(iv, iv+1);
107  targetCubPoints[0][iv] = view_type("targetCubPoints",1,dim);
108  CellTools<SpT>::getReferenceVertex(coord, cellTopo, iv);
109  for(ordinal_type d=0; d<dim; d++) {
110  basisCubPoints[0][iv](0,d) = coord(d);
111  targetCubPoints[0][iv](0,d) = coord(d);
112  }
113  }
114 
115  if ((name == "Intrepid2_HCURL_QUAD_In_FEM") || (name == "Intrepid2_HCURL_HEX_In_FEM") ||
116  (name == "Intrepid2_HCURL_TRI_In_FEM") || (name == "Intrepid2_HCURL_TET_In_FEM")) {
117  edgeBasisCubDegree = basisCubDegree - 1;
118  faceBasisCubDegree = basisCubDegree;
119  }
120  else if ((name == "Intrepid2_HDIV_QUAD_In_FEM") || (name == "Intrepid2_HDIV_HEX_In_FEM") ||
121  (name == "Intrepid2_HDIV_TRI_In_FEM") || (name == "Intrepid2_HDIV_TET_In_FEM")) {
122  edgeBasisCubDegree = basisCubDegree - 1;
123  faceBasisCubDegree = basisCubDegree -1;
124  }
125  else {
126  edgeBasisCubDegree = basisCubDegree;
127  faceBasisCubDegree = basisCubDegree;
128  }
129 
130  DefaultCubatureFactory cub_factory;
131  for(ordinal_type ie=0; ie<numEdges; ++ie) {
132  ordinal_type cub_degree = 2*edgeBasisCubDegree;
133  subCellTopologyKey(edgeDim,ie) = cellBasis->getBaseCellTopology().getKey(edgeDim, ie);
134  auto edgeBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(cellBasis->getBaseCellTopology().getKey(edgeDim, ie), cub_degree);
135  basisPointsRange(edgeDim,ie) = range_type(numBasisEvalPoints, numBasisEvalPoints+edgeBasisCub->getNumPoints());
136  numBasisEvalPoints += edgeBasisCub->getNumPoints();
137  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, edgeBasisCub->getNumPoints());
138  basisCubPoints[edgeDim][ie] = view_type("basisCubPoints",edgeBasisCub->getNumPoints(),edgeDim);
139  basisCubWeights[edgeDim][ie] = view_type("basisCubWeights",edgeBasisCub->getNumPoints());
140  edgeBasisCub->getCubature(basisCubPoints[edgeDim][ie], basisCubWeights[edgeDim][ie]);
141 
142  cub_degree = edgeBasisCubDegree + targetCubDegree;
143  auto edgeTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(cellBasis->getBaseCellTopology().getKey(edgeDim, ie), cub_degree);
144  targetPointsRange(edgeDim,ie) = range_type(numTargetEvalPoints, numTargetEvalPoints+edgeTargetCub->getNumPoints());
145  numTargetEvalPoints += edgeTargetCub->getNumPoints();
146  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, edgeTargetCub->getNumPoints());
147  targetCubPoints[edgeDim][ie] = view_type("targetCubPoints",edgeTargetCub->getNumPoints(),edgeDim);
148  targetCubWeights[edgeDim][ie] = view_type("targetCubWeights",edgeTargetCub->getNumPoints());
149  edgeTargetCub->getCubature(targetCubPoints[edgeDim][ie], targetCubWeights[edgeDim][ie]);
150  }
151 
152  for(ordinal_type iface=0; iface<numFaces; ++iface) {
153  ordinal_type cub_degree = 2*faceBasisCubDegree;
154  subCellTopologyKey(faceDim,iface) = cellBasis->getBaseCellTopology().getKey(faceDim, iface);
155  auto faceBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
156  basisPointsRange(faceDim,iface) = range_type(numBasisEvalPoints, numBasisEvalPoints+faceBasisCub->getNumPoints());
157  numBasisEvalPoints += faceBasisCub->getNumPoints();
158  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, faceBasisCub->getNumPoints());
159  basisCubPoints[faceDim][iface] = view_type("basisCubPoints",faceBasisCub->getNumPoints(),faceDim);
160  basisCubWeights[faceDim][iface] = view_type("basisCubWeights",faceBasisCub->getNumPoints());
161  faceBasisCub->getCubature(basisCubPoints[faceDim][iface], basisCubWeights[faceDim][iface]);
162 
163  cub_degree = faceBasisCubDegree + targetCubDegree;
164  auto faceTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
165  targetPointsRange(faceDim,iface) = range_type(numTargetEvalPoints, numTargetEvalPoints+faceTargetCub->getNumPoints());
166  numTargetEvalPoints += faceTargetCub->getNumPoints();
167  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, faceTargetCub->getNumPoints());
168  targetCubPoints[faceDim][iface] = view_type("targetCubPoints",faceTargetCub->getNumPoints(),faceDim);
169  targetCubWeights[faceDim][iface] = view_type("targetCubWeights",faceTargetCub->getNumPoints());
170  faceTargetCub->getCubature(targetCubPoints[faceDim][iface], targetCubWeights[faceDim][iface]);
171  }
172  subCellTopologyKey(dim,0) = cellBasis->getBaseCellTopology().getBaseKey();
173  if(hasCellDofs) {
174  ordinal_type cub_degree = 2*basisCubDegree;
175  auto elemBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
176  basisPointsRange(dim,0) = range_type(numBasisEvalPoints, numBasisEvalPoints+elemBasisCub->getNumPoints());
177  numBasisEvalPoints += elemBasisCub->getNumPoints();
178  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, elemBasisCub->getNumPoints());
179  basisCubPoints[dim][0] = view_type("basisCubPoints",elemBasisCub->getNumPoints(),dim);
180  basisCubWeights[dim][0] = view_type("basisCubWeights",elemBasisCub->getNumPoints());
181  elemBasisCub->getCubature(basisCubPoints[dim][0], basisCubWeights[dim][0]);
182 
183  cub_degree = basisCubDegree + targetCubDegree;
184  auto elemTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
185  targetPointsRange(dim,0) = range_type(numTargetEvalPoints, numTargetEvalPoints+elemTargetCub->getNumPoints());
186  numTargetEvalPoints += elemTargetCub->getNumPoints();
187  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, elemTargetCub->getNumPoints());
188  targetCubPoints[dim][0] = view_type("targetCubPoints",elemTargetCub->getNumPoints(),dim);
189  targetCubWeights[dim][0] = view_type("targetCubWeights",elemTargetCub->getNumPoints());
190  elemTargetCub->getCubature(targetCubPoints[dim][0], targetCubWeights[dim][0]);
191  }
192 }
193 
194 template<typename SpT, typename ValueType>
195 template<typename BasisPtrType>
197  const ordinal_type targetCubDegree,
198  const ordinal_type targetGradCubDegree) {
199  const auto cellTopo = cellBasis->getBaseCellTopology();
200  std::string name = cellBasis->getName();
201  ordinal_type dim = cellTopo.getDimension();
202  numBasisEvalPoints = 0;
203  numBasisDerivEvalPoints = 0;
204  numTargetEvalPoints = 0;
205  numTargetDerivEvalPoints = 0;
206  const ordinal_type edgeDim = 1;
207  const ordinal_type faceDim = 2;
208 
209  ordinal_type basisCubDegree = cellBasis->getDegree();
210  ordinal_type edgeBasisCubDegree = basisCubDegree;
211  ordinal_type faceBasisCubDegree = basisCubDegree;
212 
213  ordinal_type numVertices = (cellBasis->getDofCount(0, 0) > 0) ? cellTopo.getVertexCount() : 0;
214  ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount(): 0;
215  ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;
216 
217  INTREPID2_TEST_FOR_ABORT( (numFaces > maxSubCellsCount) || (numEdges > maxSubCellsCount),
218  ">>> ERROR (Intrepid2::ProjectionStruct:createHDivProjectionStruct, Projections do not support a cell topology with this cub cells count");
219 
220 
221  ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);
222 
223  maxNumBasisEvalPoints = numVertices; maxNumTargetEvalPoints = numVertices;
224  maxNumBasisDerivEvalPoints = 0; maxNumTargetDerivEvalPoints = 0;
225 
226  basisPointsRange = range_tag("basisPointsRange", numberSubCellDims,maxSubCellsCount);
227  targetPointsRange = range_tag("targetPointsRange", numberSubCellDims,maxSubCellsCount);
228  basisDerivPointsRange = range_tag("basisDerivPointsRange", numberSubCellDims,maxSubCellsCount);
229  targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
230  subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);
231 
232  numBasisEvalPoints += numVertices;
233  numTargetEvalPoints += numVertices;
234  view_type coord("vertex_coord", dim);
235  for(ordinal_type iv=0; iv<numVertices; ++iv) {
236  basisPointsRange(0,iv) = range_type(iv, iv+1);
237  basisCubPoints[0][iv] = view_type("basisCubPoints",1,dim);
238  targetPointsRange(0,iv) = range_type(iv, iv+1);
239  targetCubPoints[0][iv] = view_type("targetCubPoints",1,dim);
240  CellTools<SpT>::getReferenceVertex(coord, cellTopo, iv);
241  for(ordinal_type d=0; d<dim; d++) {
242  basisCubPoints[0][iv](0,d) = coord(d);
243  targetCubPoints[0][iv](0,d) = coord(d);
244  }
245  }
246 
247  DefaultCubatureFactory cub_factory;
248  for(ordinal_type ie=0; ie<numEdges; ++ie) {
249  ordinal_type cub_degree = 2*edgeBasisCubDegree;
250  subCellTopologyKey(edgeDim,ie) = cellBasis->getBaseCellTopology().getKey(edgeDim, ie);
251  auto edgeBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(cellBasis->getBaseCellTopology().getKey(edgeDim, ie), cub_degree);
252  basisDerivPointsRange(edgeDim,ie) = range_type(numBasisDerivEvalPoints, numBasisDerivEvalPoints+edgeBasisCub->getNumPoints());
253  numBasisDerivEvalPoints += edgeBasisCub->getNumPoints();
254  maxNumBasisDerivEvalPoints = std::max(maxNumBasisDerivEvalPoints, edgeBasisCub->getNumPoints());
255  basisDerivCubPoints[edgeDim][ie] = view_type("basisDerivCubPoints",edgeBasisCub->getNumPoints(),edgeDim);
256  basisDerivCubWeights[edgeDim][ie] = view_type("basisDerivCubWeights",edgeBasisCub->getNumPoints());
257  edgeBasisCub->getCubature(basisDerivCubPoints[edgeDim][ie], basisDerivCubWeights[edgeDim][ie]);
258 
259  cub_degree = edgeBasisCubDegree + targetGradCubDegree;
260  auto edgeTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(cellBasis->getBaseCellTopology().getKey(edgeDim, ie), cub_degree);
261  targetDerivPointsRange(edgeDim,ie) = range_type(numTargetDerivEvalPoints, numTargetDerivEvalPoints+edgeTargetCub->getNumPoints());
262  numTargetDerivEvalPoints += edgeTargetCub->getNumPoints();
263  maxNumTargetDerivEvalPoints = std::max(maxNumTargetDerivEvalPoints, edgeTargetCub->getNumPoints());
264  targetDerivCubPoints[edgeDim][ie] = view_type("targetDerivCubPoints",edgeTargetCub->getNumPoints(),edgeDim);
265  targetDerivCubWeights[edgeDim][ie] = view_type("targetDerivCubWeights",edgeTargetCub->getNumPoints());
266  edgeTargetCub->getCubature(targetDerivCubPoints[edgeDim][ie], targetDerivCubWeights[edgeDim][ie]);
267  }
268 
269  for(ordinal_type iface=0; iface<numFaces; ++iface) {
270  ordinal_type cub_degree = 2*faceBasisCubDegree;
271  subCellTopologyKey(faceDim,iface) = cellBasis->getBaseCellTopology().getKey(faceDim, iface);
272  auto faceBasisGradCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
273  basisDerivPointsRange(faceDim,iface) = range_type(numBasisDerivEvalPoints, numBasisDerivEvalPoints+faceBasisGradCub->getNumPoints());
274  numBasisDerivEvalPoints += faceBasisGradCub->getNumPoints();
275  maxNumBasisDerivEvalPoints = std::max(maxNumBasisDerivEvalPoints, faceBasisGradCub->getNumPoints());
276  basisDerivCubPoints[faceDim][iface] = view_type("basisDerivCubPoints",faceBasisGradCub->getNumPoints(),faceDim);
277  basisDerivCubWeights[faceDim][iface] = view_type("basisDerivCubWeights",faceBasisGradCub->getNumPoints());
278  faceBasisGradCub->getCubature(basisDerivCubPoints[faceDim][iface], basisDerivCubWeights[faceDim][iface]);
279 
280  cub_degree = faceBasisCubDegree + targetGradCubDegree;
281  auto faceTargetDerivCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
282  targetDerivPointsRange(faceDim,iface) = range_type(numTargetDerivEvalPoints, numTargetDerivEvalPoints+faceTargetDerivCub->getNumPoints());
283  numTargetDerivEvalPoints += faceTargetDerivCub->getNumPoints();
284  maxNumTargetDerivEvalPoints = std::max(maxNumTargetDerivEvalPoints, faceTargetDerivCub->getNumPoints());
285  targetDerivCubPoints[faceDim][iface] = view_type("targetDerivCubPoints",faceTargetDerivCub->getNumPoints(),faceDim);
286  targetDerivCubWeights[faceDim][iface] = view_type("targetDerivCubWeights",faceTargetDerivCub->getNumPoints());
287  faceTargetDerivCub->getCubature(targetDerivCubPoints[faceDim][iface], targetDerivCubWeights[faceDim][iface]);
288  }
289  subCellTopologyKey(dim,0) = cellBasis->getBaseCellTopology().getBaseKey();
290  if(hasCellDofs) {
291  ordinal_type cub_degree = 2*basisCubDegree;
292  auto elemBasisGradCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
293  basisDerivPointsRange(dim,0) = range_type(numBasisDerivEvalPoints, numBasisDerivEvalPoints+elemBasisGradCub->getNumPoints());
294  numBasisDerivEvalPoints += elemBasisGradCub->getNumPoints();
295  maxNumBasisDerivEvalPoints = std::max(maxNumBasisDerivEvalPoints, elemBasisGradCub->getNumPoints());
296  basisDerivCubPoints[dim][0] = view_type("basisDerivCubPoints",elemBasisGradCub->getNumPoints(),dim);
297  basisDerivCubWeights[dim][0] = view_type("basisDerivCubWeights",elemBasisGradCub->getNumPoints());
298  elemBasisGradCub->getCubature(basisDerivCubPoints[dim][0], basisDerivCubWeights[dim][0]);
299 
300  cub_degree = basisCubDegree + targetGradCubDegree;
301  auto elemTargetGradCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
302  targetDerivPointsRange(dim,0) = range_type(numTargetDerivEvalPoints, numTargetDerivEvalPoints+elemTargetGradCub->getNumPoints());
303  numTargetDerivEvalPoints += elemTargetGradCub->getNumPoints();
304  maxNumTargetDerivEvalPoints = std::max(maxNumTargetDerivEvalPoints, elemTargetGradCub->getNumPoints());
305  targetDerivCubPoints[dim][0] = view_type("targetDerivCubPoints",elemTargetGradCub->getNumPoints(),dim);
306  targetDerivCubWeights[dim][0] = view_type("targetDerivCubWeights",elemTargetGradCub->getNumPoints());
307  elemTargetGradCub->getCubature(targetDerivCubPoints[dim][0], targetDerivCubWeights[dim][0]);
308  }
309 }
310 
311 template<typename SpT, typename ValueType>
312 template<typename BasisPtrType>
314  const ordinal_type targetCubDegree,
315  const ordinal_type targetCurlCubDegre) {
316  const auto cellTopo = cellBasis->getBaseCellTopology();
317  std::string name = cellBasis->getName();
318  ordinal_type dim = cellTopo.getDimension();
319  numBasisEvalPoints = 0;
320  numBasisDerivEvalPoints = 0;
321  numTargetEvalPoints = 0;
322  numTargetDerivEvalPoints = 0;
323  const ordinal_type edgeDim = 1;
324  const ordinal_type faceDim = 2;
325 
326  ordinal_type basisCubDegree = cellBasis->getDegree();
327  ordinal_type edgeBasisCubDegree = basisCubDegree - 1;
328  ordinal_type faceBasisCubDegree = basisCubDegree;
329  ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount() : 0;
330  ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;
331  ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);
332 
333  maxNumBasisEvalPoints = 0; maxNumTargetEvalPoints = 0;
334  maxNumBasisDerivEvalPoints = 0; maxNumTargetDerivEvalPoints = 0;
335 
336  basisPointsRange = range_tag("basisPointsRange", numberSubCellDims,maxSubCellsCount);
337  targetPointsRange = range_tag("targetPointsRange", numberSubCellDims,maxSubCellsCount);
338  basisDerivPointsRange = range_tag("basisDerivPointsRange", numberSubCellDims,maxSubCellsCount);
339  targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
340  subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);
341 
342  DefaultCubatureFactory cub_factory;
343  for(ordinal_type ie=0; ie<numEdges; ++ie) {
344  ordinal_type cub_degree = 2*edgeBasisCubDegree;
345  subCellTopologyKey(edgeDim,ie) = cellBasis->getBaseCellTopology().getKey(edgeDim, ie);
346  auto edgeBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(edgeDim,ie), cub_degree);
347  basisPointsRange(edgeDim,ie) = range_type(numBasisEvalPoints, numBasisEvalPoints+edgeBasisCub->getNumPoints());
348  numBasisEvalPoints += edgeBasisCub->getNumPoints();
349  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, edgeBasisCub->getNumPoints());
350  basisCubPoints[edgeDim][ie] = view_type("basisCubPoints",edgeBasisCub->getNumPoints(),edgeDim);
351  basisCubWeights[edgeDim][ie] = view_type("basisCubWeights",edgeBasisCub->getNumPoints());
352  edgeBasisCub->getCubature(basisCubPoints[edgeDim][ie], basisCubWeights[edgeDim][ie]);
353 
354  cub_degree = edgeBasisCubDegree + targetCubDegree;
355  auto edgeTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(edgeDim,ie), cub_degree);
356  targetPointsRange(edgeDim,ie) = range_type(numTargetEvalPoints, numTargetEvalPoints+edgeTargetCub->getNumPoints());
357  numTargetEvalPoints += edgeTargetCub->getNumPoints();
358  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, edgeTargetCub->getNumPoints());
359  targetCubPoints[edgeDim][ie] = view_type("targetCubPoints",edgeTargetCub->getNumPoints(),edgeDim);
360  targetCubWeights[edgeDim][ie] = view_type("targetCubWeights",edgeTargetCub->getNumPoints());
361  edgeTargetCub->getCubature(targetCubPoints[edgeDim][ie], targetCubWeights[edgeDim][ie]);
362  }
363 
364  for(ordinal_type iface=0; iface<numFaces; ++iface) {
365  ordinal_type cub_degree = 2*faceBasisCubDegree;
366  subCellTopologyKey(faceDim,iface) = cellBasis->getBaseCellTopology().getKey(faceDim, iface);
367  auto faceBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
368  basisPointsRange(faceDim,iface) = range_type(numBasisEvalPoints, numBasisEvalPoints+faceBasisCub->getNumPoints());
369  numBasisEvalPoints += faceBasisCub->getNumPoints();
370  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, faceBasisCub->getNumPoints());
371  basisCubPoints[faceDim][iface] = view_type("basisCubPoints",faceBasisCub->getNumPoints(),faceDim);
372  basisCubWeights[faceDim][iface] = view_type("basisCubWeights",faceBasisCub->getNumPoints());
373  faceBasisCub->getCubature(basisCubPoints[faceDim][iface], basisCubWeights[faceDim][iface]);
374 
375  auto faceBasisDerivCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
376  basisDerivPointsRange(faceDim,iface) = range_type(numBasisDerivEvalPoints, numBasisDerivEvalPoints+faceBasisCub->getNumPoints());
377  numBasisDerivEvalPoints += faceBasisCub->getNumPoints();
378  maxNumBasisDerivEvalPoints = std::max(maxNumBasisDerivEvalPoints, faceBasisCub->getNumPoints());
379  basisDerivCubPoints[faceDim][iface] = view_type("basisDerivCubPoints",faceBasisCub->getNumPoints(),faceDim);
380  basisDerivCubWeights[faceDim][iface] = view_type("basisDerivCubWeights",faceBasisCub->getNumPoints());
381  faceBasisCub->getCubature(basisDerivCubPoints[faceDim][iface], basisDerivCubWeights[faceDim][iface]);
382 
383  cub_degree = faceBasisCubDegree + targetCubDegree;
384  auto faceTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
385  targetPointsRange(faceDim,iface) = range_type(numTargetEvalPoints, numTargetEvalPoints+faceTargetCub->getNumPoints());
386  numTargetEvalPoints += faceTargetCub->getNumPoints();
387  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, faceTargetCub->getNumPoints());
388  targetCubPoints[faceDim][iface] = view_type("targetCubPoints",faceTargetCub->getNumPoints(),faceDim);
389  targetCubWeights[faceDim][iface] = view_type("targetCubWeights",faceTargetCub->getNumPoints());
390  faceTargetCub->getCubature(targetCubPoints[faceDim][iface], targetCubWeights[faceDim][iface]);
391 
392  cub_degree = faceBasisCubDegree + targetCurlCubDegre;
393  auto faceTargetDerivCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(faceDim,iface), cub_degree);
394  targetDerivPointsRange(faceDim,iface) = range_type(numTargetDerivEvalPoints, numTargetDerivEvalPoints+faceTargetDerivCub->getNumPoints());
395  numTargetDerivEvalPoints += faceTargetDerivCub->getNumPoints();
396  maxNumTargetDerivEvalPoints = std::max(maxNumTargetDerivEvalPoints, faceTargetDerivCub->getNumPoints());
397  targetDerivCubPoints[faceDim][iface] = view_type("targetDerivCubPoints",faceTargetDerivCub->getNumPoints(),faceDim);
398  targetDerivCubWeights[faceDim][iface] = view_type("targetDerivCubWeights",faceTargetDerivCub->getNumPoints());
399  faceTargetDerivCub->getCubature(targetDerivCubPoints[faceDim][iface], targetDerivCubWeights[faceDim][iface]);
400  }
401 
402  subCellTopologyKey(dim,0) = cellBasis->getBaseCellTopology().getBaseKey();
403  if(hasCellDofs) {
404  ordinal_type cub_degree = 2*basisCubDegree;
405  auto elemBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
406  basisPointsRange(dim,0) = range_type(numBasisEvalPoints, numBasisEvalPoints+elemBasisCub->getNumPoints());
407  numBasisEvalPoints += elemBasisCub->getNumPoints();
408  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, elemBasisCub->getNumPoints());
409  basisCubPoints[dim][0] = view_type("basisCubPoints",elemBasisCub->getNumPoints(),dim);
410  basisCubWeights[dim][0] = view_type("basisCubWeights",elemBasisCub->getNumPoints());
411  elemBasisCub->getCubature(basisCubPoints[dim][0], basisCubWeights[dim][0]);
412 
413  basisDerivPointsRange(dim,0) = range_type(numBasisDerivEvalPoints, numBasisDerivEvalPoints+elemBasisCub->getNumPoints());
414  numBasisDerivEvalPoints += elemBasisCub->getNumPoints();
415  maxNumBasisDerivEvalPoints = std::max(maxNumBasisDerivEvalPoints, elemBasisCub->getNumPoints());
416  basisDerivCubPoints[dim][0] = view_type("basisDerivCubPoints",elemBasisCub->getNumPoints(),dim);
417  basisDerivCubWeights[dim][0] = view_type("basisDerivCubWeights",elemBasisCub->getNumPoints());
418  elemBasisCub->getCubature(basisDerivCubPoints[dim][0], basisDerivCubWeights[dim][0]);
419 
420  cub_degree = basisCubDegree + targetCubDegree;
421  auto elemTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
422  targetPointsRange(dim,0) = range_type(numTargetEvalPoints, numTargetEvalPoints+elemTargetCub->getNumPoints());
423  numTargetEvalPoints += elemTargetCub->getNumPoints();
424  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, elemTargetCub->getNumPoints());
425  targetCubPoints[dim][0] = view_type("targetCubPoints",elemTargetCub->getNumPoints(),dim);
426  targetCubWeights[dim][0] = view_type("targetCubWeights",elemTargetCub->getNumPoints());
427  elemTargetCub->getCubature(targetCubPoints[dim][0], targetCubWeights[dim][0]);
428 
429  cub_degree = basisCubDegree + targetCurlCubDegre;
430  auto elemTargetCurlCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
431  targetDerivPointsRange(dim,0) = range_type(numTargetDerivEvalPoints, numTargetDerivEvalPoints+elemTargetCurlCub->getNumPoints());
432  numTargetDerivEvalPoints += elemTargetCurlCub->getNumPoints();
433  maxNumTargetDerivEvalPoints = std::max(maxNumTargetDerivEvalPoints, elemTargetCurlCub->getNumPoints());
434  targetDerivCubPoints[dim][0] = view_type("targetDerivCubPoints",elemTargetCurlCub->getNumPoints(),dim);
435  targetDerivCubWeights[dim][0] = view_type("targetDerivCubWeights",elemTargetCurlCub->getNumPoints());
436  elemTargetCurlCub->getCubature(targetDerivCubPoints[dim][0], targetDerivCubWeights[dim][0]);
437  }
438 }
439 
440 template<typename SpT, typename ValueType>
441 template<typename BasisPtrType>
443  const ordinal_type targetCubDegree,
444  const ordinal_type targetDivCubDegre) {
445  const auto cellTopo = cellBasis->getBaseCellTopology();
446  std::string name = cellBasis->getName();
447  ordinal_type dim = cellTopo.getDimension();
448  numBasisEvalPoints = 0;
449  numBasisDerivEvalPoints = 0;
450  numTargetEvalPoints = 0;
451  numTargetDerivEvalPoints = 0;
452  const ordinal_type sideDim = dim - 1;
453 
454  ordinal_type basisCubDegree = cellBasis->getDegree();
455  ordinal_type sideBasisCubDegree = basisCubDegree - 1;
456  ordinal_type numSides = cellTopo.getSideCount()*ordinal_type(cellBasis->getDofCount(sideDim, 0) > 0);
457  ordinal_type hasCellDofs = (cellBasis->getDofCount(dim, 0) > 0);
458 
459  INTREPID2_TEST_FOR_ABORT( numSides > maxSubCellsCount,
460  ">>> ERROR (Intrepid2::ProjectionStruct:createHDivProjectionStruct, Projections do not support a cell topology with so many sides");
461 
462  maxNumBasisEvalPoints = 0; maxNumTargetEvalPoints = 0;
463  maxNumBasisDerivEvalPoints = 0; maxNumTargetDerivEvalPoints = 0;
464 
465  basisPointsRange = range_tag("basisPointsRange", numberSubCellDims,maxSubCellsCount);
466  targetPointsRange = range_tag("targetPointsRange", numberSubCellDims,maxSubCellsCount);
467  basisDerivPointsRange = range_tag("basisDerivPointsRange", numberSubCellDims,maxSubCellsCount);
468  targetDerivPointsRange = range_tag("targetDerivPointsRange", numberSubCellDims,maxSubCellsCount);
469  subCellTopologyKey = key_tag("subCellTopologyKey",numberSubCellDims,maxSubCellsCount);
470 
472  if(cellTopo.getKey() == shards::getCellTopologyData<shards::Hexahedron<8> >()->key)
473  hcurlBasis = new Basis_HCURL_HEX_In_FEM<host_space_type,ValueType,ValueType>(cellBasis->getDegree());
474  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Tetrahedron<4> >()->key)
475  hcurlBasis = new Basis_HCURL_TET_In_FEM<host_space_type,ValueType,ValueType>(cellBasis->getDegree());
476  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Quadrilateral<4> >()->key)
477  hcurlBasis = new Basis_HGRAD_QUAD_Cn_FEM<host_space_type,ValueType,ValueType>(cellBasis->getDegree());
478  else if(cellTopo.getKey() == shards::getCellTopologyData<shards::Triangle<3> >()->key)
479  hcurlBasis = new Basis_HGRAD_TRI_Cn_FEM<host_space_type,ValueType,ValueType>(cellBasis->getDegree());
480  else {
481  std::stringstream ss;
482  ss << ">>> ERROR (Intrepid2::ProjectionTools::createHDivProjectionStruct): "
483  << "Method not implemented for basis " << name;
484  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
485  }
486 
487  bool haveHCurlConstraint = (hcurlBasis != NULL) && (hcurlBasis->getDofCount(dim,0) > 0);
488  if(hcurlBasis != NULL) delete hcurlBasis;
489 
490  DefaultCubatureFactory cub_factory;
491 
492  for(ordinal_type is=0; is<numSides; ++is) {
493  ordinal_type cub_degree = 2*sideBasisCubDegree;
494  subCellTopologyKey(sideDim,is) = cellBasis->getBaseCellTopology().getKey(sideDim, is);
495  auto sideBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(sideDim,is), cub_degree);
496  basisPointsRange(sideDim,is) = range_type(numBasisEvalPoints, numBasisEvalPoints+sideBasisCub->getNumPoints());
497  numBasisEvalPoints += sideBasisCub->getNumPoints();
498  basisCubPoints[sideDim][is] = view_type("basisCubPoints",sideBasisCub->getNumPoints(),sideDim);
499  basisCubWeights[sideDim][is] = view_type("basisCubWeights",sideBasisCub->getNumPoints());
500  sideBasisCub->getCubature(basisCubPoints[sideDim][is], basisCubWeights[sideDim][is]);
501  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, sideBasisCub->getNumPoints());
502 
503  cub_degree = sideBasisCubDegree + targetCubDegree;
504  auto sideTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(sideDim,is), cub_degree);
505  targetPointsRange(sideDim,is) = range_type(numTargetEvalPoints, numTargetEvalPoints+sideTargetCub->getNumPoints());
506  numTargetEvalPoints += sideTargetCub->getNumPoints();
507  targetCubPoints[sideDim][is] = view_type("targetCubPoints",sideTargetCub->getNumPoints(),sideDim);
508  targetCubWeights[sideDim][is] = view_type("targetCubWeights",sideTargetCub->getNumPoints());
509  sideTargetCub->getCubature(targetCubPoints[sideDim][is], targetCubWeights[sideDim][is]);
510  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, sideTargetCub->getNumPoints());
511  }
512 
513  subCellTopologyKey(dim,0) = cellBasis->getBaseCellTopology().getBaseKey();
514  if(hasCellDofs) {
515  ordinal_type cub_degree = 2*basisCubDegree - 1;
516  auto elemBasisDivCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
517  basisDerivPointsRange(dim,0) = range_type(numBasisDerivEvalPoints, numBasisDerivEvalPoints+elemBasisDivCub->getNumPoints());
518  numBasisDerivEvalPoints += elemBasisDivCub->getNumPoints();
519  basisDerivCubPoints[dim][0] = view_type("basisDerivCubPoints",elemBasisDivCub->getNumPoints(),dim);
520  basisDerivCubWeights[dim][0] = view_type("basisDerivCubWeights",elemBasisDivCub->getNumPoints());
521  elemBasisDivCub->getCubature(basisDerivCubPoints[dim][0], basisDerivCubWeights[dim][0]);
522  maxNumBasisDerivEvalPoints = std::max(maxNumBasisDerivEvalPoints, elemBasisDivCub->getNumPoints());
523 
524  cub_degree = basisCubDegree - 1 + targetDivCubDegre;
525  auto elemTargetDivCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
526  targetDerivPointsRange(dim,0) = range_type(numTargetDerivEvalPoints, numTargetDerivEvalPoints+elemTargetDivCub->getNumPoints());
527  numTargetDerivEvalPoints += elemTargetDivCub->getNumPoints();
528  targetDerivCubPoints[dim][0] = view_type("targetDerivCubPoints",elemTargetDivCub->getNumPoints(),dim);
529  targetDerivCubWeights[dim][0] = view_type("targetDerivCubWeights",elemTargetDivCub->getNumPoints());
530  elemTargetDivCub->getCubature(targetDerivCubPoints[dim][0], targetDerivCubWeights[dim][0]);
531  maxNumTargetDerivEvalPoints = std::max(maxNumTargetDerivEvalPoints, elemTargetDivCub->getNumPoints());
532 
533  if(haveHCurlConstraint)
534  {
535  cub_degree = 2*basisCubDegree;
536  auto elemBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
537  basisPointsRange(dim,0) = range_type(numBasisEvalPoints, numBasisEvalPoints + elemBasisCub->getNumPoints());
538  numBasisEvalPoints += elemBasisCub->getNumPoints();
539  basisCubPoints[dim][0] = view_type("basisCubPoints",elemBasisCub->getNumPoints(),dim);
540  basisCubWeights[dim][0] = view_type("basisCubWeights",elemBasisCub->getNumPoints());
541  elemBasisCub->getCubature(basisCubPoints[dim][0], basisCubWeights[dim][0]);
542  maxNumBasisEvalPoints = std::max(maxNumBasisEvalPoints, elemBasisCub->getNumPoints());
543 
544  cub_degree = basisCubDegree + targetCubDegree;
545  auto elemTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(subCellTopologyKey(dim,0), cub_degree);
546  targetPointsRange(dim,0) = range_type(numTargetEvalPoints, numTargetEvalPoints + elemTargetCub->getNumPoints());
547  numTargetEvalPoints += elemTargetCub->getNumPoints();
548  targetCubPoints[dim][0] = view_type("targetCubPoints",elemTargetCub->getNumPoints(),dim);
549  targetCubWeights[dim][0] = view_type("targetCubWeights",elemTargetCub->getNumPoints());
550  elemTargetCub->getCubature(targetCubPoints[dim][0], targetCubWeights[dim][0]);
551  maxNumTargetEvalPoints = std::max(maxNumTargetEvalPoints, elemTargetCub->getNumPoints());
552  }
553  }
554 }
555 
556 template<typename SpT, typename ValueType>
557 template<typename BasisPtrType>
559  const ordinal_type targetCubDegree) {
560  const auto cellTopo = cellBasis->getBaseCellTopology();
561  ordinal_type dim = cellTopo.getDimension();
562  numBasisEvalPoints = 0;
563  numBasisDerivEvalPoints = 0;
564  numTargetEvalPoints = 0;
565  numTargetDerivEvalPoints = 0;
566 
567  basisPointsRange = range_tag("basisPointsRange", 4,maxSubCellsCount);
568  targetPointsRange = range_tag("targetPointsRange", 4,maxSubCellsCount);
569  basisDerivPointsRange = range_tag("basisDerivPointsRange", 4,maxSubCellsCount);
570  targetDerivPointsRange = range_tag("targetDerivPointsRange", 4,maxSubCellsCount);
571  subCellTopologyKey = key_tag("subCellTopologyKey",4,maxSubCellsCount);
572 
573  ordinal_type basisCubDegree = cellBasis->getDegree();
574  DefaultCubatureFactory cub_factory;
575  subCellTopologyKey(dim,0) = cellBasis->getBaseCellTopology().getBaseKey();
576 
577  maxNumBasisEvalPoints = 0; maxNumTargetEvalPoints =0;
578  if(cellBasis->getDofCount(dim,0)>0) {
579 
580  ordinal_type cub_degree = 2*basisCubDegree;
581  auto elemBasisCub = cub_factory.create<host_space_type, ValueType, ValueType>(cellTopo.getBaseKey(), cub_degree);
582  basisPointsRange(dim,0) = range_type(0, elemBasisCub->getNumPoints());
583  numBasisEvalPoints += elemBasisCub->getNumPoints();
584  maxNumBasisEvalPoints = elemBasisCub->getNumPoints();
585  basisCubPoints[dim][0] = view_type("basisCubPoints",elemBasisCub->getNumPoints(),dim);
586  basisCubWeights[dim][0] = view_type("basisCubWeights",elemBasisCub->getNumPoints());
587  elemBasisCub->getCubature(basisCubPoints[dim][0], basisCubWeights[dim][0]);
588 
589  cub_degree = basisCubDegree + targetCubDegree;
590  auto elemTargetCub = cub_factory.create<host_space_type, ValueType, ValueType>(cellTopo.getBaseKey(), cub_degree);
591  targetPointsRange(dim,0) = range_type(0, elemTargetCub->getNumPoints());
592  numTargetEvalPoints += elemTargetCub->getNumPoints();
593  maxNumTargetEvalPoints = elemTargetCub->getNumPoints();
594  targetCubPoints[dim][0] = view_type("targetCubPoints",elemTargetCub->getNumPoints(),dim);
595  targetCubWeights[dim][0] = view_type("targetCubWeights",elemTargetCub->getNumPoints());
596  elemTargetCub->getCubature(targetCubPoints[dim][0], targetCubWeights[dim][0]);
597  }
598 }
599 
600 }
601 }
602 #endif
603 
604 
605 
606 
607 
Header file for the abstract base class Intrepid2::DefaultCubatureFactory.
Header file for the Intrepid2::Basis_HCURL_HEX_In_FEM class.
A factory class that generates specific instances of cubatures.
void createL2ProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
Initialize the ProjectionStruct for L2 projections.
void createHDivProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetDivCubDegre)
Initialize the ProjectionStruct for HDIV projections.
Implementation of the default H(curl)-compatible Nedelec (first kind) basis of arbitrary degree on Te...
static void getReferenceVertex(Kokkos::DynRankView< cellVertexValueType, cellVertexProperties...> cellVertex, const shards::CellTopology cell, const ordinal_type vertexOrd)
Retrieves the Cartesian coordinates of a reference cell vertex.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
Implementation of the default H(curl)-compatible FEM basis on Hexahedron cell.
Header function for Intrepid2::Util class and other utility functions.
static Teuchos::RCP< Cubature< ExecSpaceType, pointValueType, weightValueType > > create(unsigned topologyKey, const std::vector< ordinal_type > &degree, const EPolyType polytype=POLYTYPE_MAX)
Factory method.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_In_FEM class.
ordinal_type getDofCount(const ordinal_type subcDim, const ordinal_type subcOrd) const
DoF count for specified subcell.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
void createHVolProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
Initialize the ProjectionStruct for HVOL (local-L2) projection.
void createHCurlProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetCurlCubDegre)
Initialize the ProjectionStruct for HCURL projections.
void createHGradProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetGradCubDegre)
Initialize the ProjectionStruct for HGRAD projections.
Implementation of the default H(grad)-compatible Lagrange basis of arbitrary degree on Triangle cell...
Implementation of the default H(grad)-compatible FEM basis of degree n on Quadrilateral cell Implemen...