92struct ComputeBasisCoeffsOnEdges_L2 {
93 const ViewType1 basisCoeffs_;
94 const ViewType2 negPartialProj_;
95 const ViewType2 basisDofDofAtBasisEPoints_;
96 const ViewType2 basisAtBasisEPoints_;
97 const ViewType3 basisEWeights_;
98 const ViewType2 wBasisDofAtBasisEPoints_;
99 const ViewType3 targetEWeights_;
100 const ViewType2 basisAtTargetEPoints_;
101 const ViewType2 wBasisDofAtTargetEPoints_;
102 const ViewType4 computedDofs_;
103 const ViewType5 tagToOrdinal_;
104 const ViewType6 targetAtTargetEPoints_;
105 const ViewType2 targetTanAtTargetEPoints_;
106 const ViewType2 refEdgesVec_;
107 ordinal_type fieldDim_;
108 ordinal_type edgeCardinality_;
109 ordinal_type offsetBasis_;
110 ordinal_type offsetTarget_;
111 ordinal_type numVertexDofs_;
112 ordinal_type edgeDim_;
115 ComputeBasisCoeffsOnEdges_L2(
const ViewType1 basisCoeffs, ViewType2 negPartialProj,
const ViewType2 basisDofDofAtBasisEPoints,
116 const ViewType2 basisAtBasisEPoints,
const ViewType3 basisEWeights,
const ViewType2 wBasisDofAtBasisEPoints,
const ViewType3 targetEWeights,
117 const ViewType2 basisAtTargetEPoints,
const ViewType2 wBasisDofAtTargetEPoints,
const ViewType4 computedDofs,
const ViewType5 tagToOrdinal,
118 const ViewType6 targetAtTargetEPoints,
const ViewType2 targetTanAtTargetEPoints,
const ViewType2 refEdgesVec,
119 ordinal_type fieldDim, ordinal_type edgeCardinality, ordinal_type offsetBasis,
120 ordinal_type offsetTarget, ordinal_type numVertexDofs, ordinal_type edgeDim, ordinal_type iedge) :
121 basisCoeffs_(basisCoeffs), negPartialProj_(negPartialProj), basisDofDofAtBasisEPoints_(basisDofDofAtBasisEPoints),
122 basisAtBasisEPoints_(basisAtBasisEPoints), basisEWeights_(basisEWeights), wBasisDofAtBasisEPoints_(wBasisDofAtBasisEPoints), targetEWeights_(targetEWeights),
123 basisAtTargetEPoints_(basisAtTargetEPoints), wBasisDofAtTargetEPoints_(wBasisDofAtTargetEPoints),
124 computedDofs_(computedDofs), tagToOrdinal_(tagToOrdinal), targetAtTargetEPoints_(targetAtTargetEPoints),
125 targetTanAtTargetEPoints_(targetTanAtTargetEPoints), refEdgesVec_(refEdgesVec),
126 fieldDim_(fieldDim), edgeCardinality_(edgeCardinality), offsetBasis_(offsetBasis),
127 offsetTarget_(offsetTarget), numVertexDofs_(numVertexDofs), edgeDim_(edgeDim), iedge_(iedge)
131 KOKKOS_INLINE_FUNCTION
132 operator()(
const ordinal_type ic)
const {
133 for(ordinal_type j=0; j <edgeCardinality_; ++j) {
134 ordinal_type jdof =tagToOrdinal_(edgeDim_, iedge_, j);
135 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq) {
136 for(ordinal_type d=0; d <fieldDim_; ++d)
137 basisDofDofAtBasisEPoints_(ic,j,iq) += basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d)*refEdgesVec_(iedge_,d);
138 wBasisDofAtBasisEPoints_(ic,j,iq) = basisDofDofAtBasisEPoints_(ic,j,iq)*basisEWeights_(iq);
140 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq) {
141 for(ordinal_type d=0; d <fieldDim_; ++d)
142 wBasisDofAtTargetEPoints_(ic,j,iq) += basisAtTargetEPoints_(ic,jdof,offsetTarget_+iq,d)*refEdgesVec_(iedge_,d)*targetEWeights_(iq);
146 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq)
147 for(ordinal_type d=0; d <fieldDim_; ++d)
148 targetTanAtTargetEPoints_(ic,iq) += targetAtTargetEPoints_(ic,offsetTarget_+iq,d)*refEdgesVec_(iedge_,d);
150 for(ordinal_type j=0; j <numVertexDofs_; ++j) {
151 ordinal_type jdof = computedDofs_(j);
152 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq)
153 for(ordinal_type d=0; d <fieldDim_; ++d)
154 negPartialProj_(ic,iq) -= basisCoeffs_(ic,jdof)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d)*refEdgesVec_(iedge_,d);
161struct ComputeBasisCoeffsOnFaces_L2 {
162 const ViewType1 basisCoeffs_;
163 const ViewType2 negPartialProj_;
164 const ViewType2 faceBasisDofAtBasisEPoints_;
165 const ViewType2 basisAtBasisEPoints_;
166 const ViewType3 basisEWeights_;
167 const ViewType2 wBasisDofAtBasisEPoints_;
168 const ViewType3 targetEWeights_;
169 const ViewType2 basisAtTargetEPoints_;
170 const ViewType2 wBasisDofAtTargetEPoints_;
171 const ViewType4 computedDofs_;
172 const ViewType5 tagToOrdinal_;
173 const ViewType6 orts_;
174 const ViewType7 targetAtTargetEPoints_;
175 const ViewType2 targetDofAtTargetEPoints_;
176 const ViewType8 faceParametrization_;
177 ordinal_type fieldDim_;
178 ordinal_type faceCardinality_;
179 ordinal_type offsetBasis_;
180 ordinal_type offsetTarget_;
181 ordinal_type numVertexEdgeDofs_;
182 ordinal_type numFaces_;
183 ordinal_type faceDim_;
187 bool isHCurlBasis_, isHDivBasis_;
189 ComputeBasisCoeffsOnFaces_L2(
const ViewType1 basisCoeffs, ViewType2 negPartialProj,
const ViewType2 faceBasisDofAtBasisEPoints,
190 const ViewType2 basisAtBasisEPoints,
const ViewType3 basisEWeights,
const ViewType2 wBasisDofAtBasisEPoints,
const ViewType3 targetEWeights,
191 const ViewType2 basisAtTargetEPoints,
const ViewType2 wBasisDofAtTargetEPoints,
const ViewType4 computedDofs,
const ViewType5 tagToOrdinal,
192 const ViewType6 orts,
const ViewType7 targetAtTargetEPoints,
const ViewType2 targetDofAtTargetEPoints,
193 const ViewType8 faceParametrization, ordinal_type fieldDim, ordinal_type faceCardinality, ordinal_type offsetBasis,
194 ordinal_type offsetTarget, ordinal_type numVertexEdgeDofs, ordinal_type numFaces, ordinal_type faceDim,
195 ordinal_type dim, ordinal_type iface,
unsigned topoKey,
bool isHCurlBasis,
bool isHDivBasis) :
196 basisCoeffs_(basisCoeffs), negPartialProj_(negPartialProj), faceBasisDofAtBasisEPoints_(faceBasisDofAtBasisEPoints),
197 basisAtBasisEPoints_(basisAtBasisEPoints), basisEWeights_(basisEWeights), wBasisDofAtBasisEPoints_(wBasisDofAtBasisEPoints), targetEWeights_(targetEWeights),
198 basisAtTargetEPoints_(basisAtTargetEPoints), wBasisDofAtTargetEPoints_(wBasisDofAtTargetEPoints),
199 computedDofs_(computedDofs), tagToOrdinal_(tagToOrdinal), orts_(orts), targetAtTargetEPoints_(targetAtTargetEPoints),
200 targetDofAtTargetEPoints_(targetDofAtTargetEPoints), faceParametrization_(faceParametrization),
201 fieldDim_(fieldDim), faceCardinality_(faceCardinality), offsetBasis_(offsetBasis),
202 offsetTarget_(offsetTarget), numVertexEdgeDofs_(numVertexEdgeDofs), numFaces_(numFaces),
203 faceDim_(faceDim), dim_(dim), iface_(iface), topoKey_(topoKey),
204 isHCurlBasis_(isHCurlBasis), isHDivBasis_(isHDivBasis)
208 KOKKOS_INLINE_FUNCTION
209 operator()(
const ordinal_type ic)
const {
211 ordinal_type fOrt[6];
212 orts_(ic).getFaceOrientation(fOrt, numFaces_);
213 ordinal_type ort = fOrt[iface_];
216 typename ViewType3::value_type data[2*3];
217 auto tangents = ViewType3(data, 2, dim_);
219 typename ViewType3::value_type tmp[2] = {};
220 for(ordinal_type j=0; j <faceCardinality_; ++j) {
221 ordinal_type jdof = tagToOrdinal_(faceDim_, iface_, j);
222 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq) {
224 for(ordinal_type d=0; d <fieldDim_; ++d) {
225 tmp[0] += tangents(0,d)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
226 tmp[1] += tangents(1,d)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
229 for(ordinal_type d=0; d <fieldDim_; ++d) {
230 faceBasisDofAtBasisEPoints_(ic,j,iq,d) = tangents(0,d)*tmp[1]-tangents(1,d)*tmp[0];
231 wBasisDofAtBasisEPoints_(ic,j,iq,d) = faceBasisDofAtBasisEPoints_(ic,j,iq,d) * basisEWeights_(iq);
234 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq) {
236 for(ordinal_type d=0; d <fieldDim_; ++d) {
237 tmp[0] += tangents(0,d)*basisAtTargetEPoints_(ic,jdof,offsetTarget_+iq,d);
238 tmp[1] += tangents(1,d)*basisAtTargetEPoints_(ic,jdof,offsetTarget_+iq,d);
241 for(ordinal_type d=0; d <fieldDim_; ++d)
242 wBasisDofAtTargetEPoints_(ic,j,iq,d) = (tangents(0,d)*tmp[1]-tangents(1,d)*tmp[0]) * targetEWeights_(iq);
246 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq) {
248 for(ordinal_type d=0; d <fieldDim_; ++d) {
249 tmp[0] += tangents(0,d)*targetAtTargetEPoints_(ic,offsetTarget_+iq,d);
250 tmp[1] += tangents(1,d)*targetAtTargetEPoints_(ic,offsetTarget_+iq,d);
253 for(ordinal_type d=0; d <fieldDim_; ++d)
254 targetDofAtTargetEPoints_(ic,iq,d) = (tangents(0,d)*tmp[1]-tangents(1,d)*tmp[0]);
257 for(ordinal_type j=0; j <numVertexEdgeDofs_; ++j) {
258 ordinal_type jdof = computedDofs_(j);
259 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq) {
261 for(ordinal_type d=0; d <fieldDim_; ++d) {
262 tmp[0] += tangents(0,d)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
263 tmp[1] += tangents(1,d)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
267 for(ordinal_type d=0; d <fieldDim_; ++d)
268 negPartialProj_(ic,iq,d) -= (tangents(0,d)*tmp[1]-tangents(1,d)*tmp[0])*basisCoeffs_(ic,jdof);
272 typename ViewType3::value_type coeff[3];
274 typename ViewType3::value_type data[3*3];
275 auto tangentsAndNormal = ViewType3(data, dim_, dim_);
277 for(ordinal_type d=0; d <dim_; ++d)
278 coeff[d] = tangentsAndNormal(dim_-1,d);
283 for(ordinal_type j=0; j <faceCardinality_; ++j) {
284 ordinal_type jdof = tagToOrdinal_(faceDim_, iface_, j);
285 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq) {
286 faceBasisDofAtBasisEPoints_(ic,j,iq,0) = 0;
287 for(ordinal_type d=0; d <fieldDim_; ++d)
288 faceBasisDofAtBasisEPoints_(ic,j,iq,0) += coeff[d]*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
289 wBasisDofAtBasisEPoints_(ic,j,iq,0) = faceBasisDofAtBasisEPoints_(ic,j,iq,0) * basisEWeights_(iq);
291 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq) {
292 typename ViewType2::value_type sum=0;
293 for(ordinal_type d=0; d <fieldDim_; ++d)
294 sum += coeff[d]*basisAtTargetEPoints_(ic,jdof,offsetTarget_+iq,d);
295 wBasisDofAtTargetEPoints_(ic,j,iq,0) = sum * targetEWeights_(iq);
299 for(ordinal_type d=0; d <fieldDim_; ++d)
300 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq)
301 targetDofAtTargetEPoints_(ic,iq,0) += coeff[d]*targetAtTargetEPoints_(ic,offsetTarget_+iq,d);
303 for(ordinal_type j=0; j <numVertexEdgeDofs_; ++j) {
304 ordinal_type jdof = computedDofs_(j);
305 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq)
306 for(ordinal_type d=0; d <fieldDim_; ++d)
307 negPartialProj_(ic,iq,0) -= basisCoeffs_(ic,jdof)*coeff[d]*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
316struct ComputeBasisCoeffsOnCells_L2 {
317 const ViewType1 basisCoeffs_;
318 const ViewType2 negPartialProj_;
319 const ViewType2 internalBasisAtBasisEPoints_;
320 const ViewType2 basisAtBasisEPoints_;
321 const ViewType3 basisEWeights_;
322 const ViewType2 wBasisAtBasisEPoints_;
323 const ViewType3 targetEWeights_;
324 const ViewType2 basisAtTargetEPoints_;
325 const ViewType2 wBasisDofAtTargetEPoints_;
326 const ViewType4 computedDofs_;
327 const ViewType5 elemDof_;
328 ordinal_type fieldDim_;
329 ordinal_type numElemDofs_;
330 ordinal_type offsetBasis_;
331 ordinal_type offsetTarget_;
332 ordinal_type numVertexEdgeFaceDofs_;
334 ComputeBasisCoeffsOnCells_L2(
const ViewType1 basisCoeffs, ViewType2 negPartialProj,
const ViewType2 internalBasisAtBasisEPoints,
335 const ViewType2 basisAtBasisEPoints,
const ViewType3 basisEWeights,
const ViewType2 wBasisAtBasisEPoints,
const ViewType3 targetEWeights,
336 const ViewType2 basisAtTargetEPoints,
const ViewType2 wBasisDofAtTargetEPoints,
const ViewType4 computedDofs,
const ViewType5 elemDof,
337 ordinal_type fieldDim, ordinal_type numElemDofs, ordinal_type offsetBasis, ordinal_type offsetTarget, ordinal_type numVertexEdgeFaceDofs) :
338 basisCoeffs_(basisCoeffs), negPartialProj_(negPartialProj), internalBasisAtBasisEPoints_(internalBasisAtBasisEPoints),
339 basisAtBasisEPoints_(basisAtBasisEPoints), basisEWeights_(basisEWeights), wBasisAtBasisEPoints_(wBasisAtBasisEPoints), targetEWeights_(targetEWeights),
340 basisAtTargetEPoints_(basisAtTargetEPoints), wBasisDofAtTargetEPoints_(wBasisDofAtTargetEPoints),
341 computedDofs_(computedDofs), elemDof_(elemDof), fieldDim_(fieldDim), numElemDofs_(numElemDofs), offsetBasis_(offsetBasis),
342 offsetTarget_(offsetTarget), numVertexEdgeFaceDofs_(numVertexEdgeFaceDofs) {}
345 KOKKOS_INLINE_FUNCTION
346 operator()(
const ordinal_type ic)
const {
348 for(ordinal_type j=0; j <numElemDofs_; ++j) {
349 ordinal_type idof = elemDof_(j);
350 for(ordinal_type d=0; d <fieldDim_; ++d) {
351 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq) {
352 internalBasisAtBasisEPoints_(ic,j,iq,d) = basisAtBasisEPoints_(ic,idof,offsetBasis_+iq,d);
353 wBasisAtBasisEPoints_(ic,j,iq,d) = internalBasisAtBasisEPoints_(ic,j,iq,d) * basisEWeights_(iq);
355 for(ordinal_type iq=0; iq <ordinal_type(targetEWeights_.extent(0)); ++iq) {
356 wBasisDofAtTargetEPoints_(ic,j,iq,d) = basisAtTargetEPoints_(ic,idof,offsetTarget_+iq,d)* targetEWeights_(iq);
360 for(ordinal_type j=0; j < numVertexEdgeFaceDofs_; ++j) {
361 ordinal_type jdof = computedDofs_(j);
362 for(ordinal_type iq=0; iq <ordinal_type(basisEWeights_.extent(0)); ++iq)
363 for(ordinal_type d=0; d <fieldDim_; ++d) {
364 negPartialProj_(ic,iq,d) -= basisCoeffs_(ic,jdof)*basisAtBasisEPoints_(ic,jdof,offsetBasis_+iq,d);
376 const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
377 const BasisType* cellBasis,
379 const EvalPointsType ePointType) {
380 const auto cellTopo = cellBasis->getBaseCellTopology();
382 ordinal_type dim = cellTopo.getDimension();
383 ordinal_type numCells = ePoints.extent(0);
384 const ordinal_type edgeDim = 1;
385 const ordinal_type faceDim = 2;
387 ordinal_type numVertices = (cellBasis->getDofCount(0, 0) > 0) ? cellTopo.getVertexCount() : 0;
388 ordinal_type numEdges = (cellBasis->getDofCount(edgeDim, 0) > 0) ? cellTopo.getEdgeCount() : 0;
389 ordinal_type numFaces = (cellBasis->getDofCount(faceDim, 0) > 0) ? cellTopo.getFaceCount() : 0;
390 ordinal_type numVols = (cellBasis->getDofCount(dim, 0) > 0);
394 typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParamEdge, subcellParamFace;
403 for(ordinal_type iv=0; iv<numVertices; ++iv) {
404 auto vertexEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getEvalPoints(0,iv,ePointType));
406 ePointsRange(0, iv), Kokkos::ALL()), vertexEPoints);
410 for(ordinal_type ie=0; ie<numEdges; ++ie) {
411 auto edgePointsRange = ePointsRange(edgeDim, ie);
412 auto edgeEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getEvalPoints(edgeDim,ie,ePointType));
414 const auto topoKey = refTopologyKey(edgeDim,ie);
417 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
418 KOKKOS_LAMBDA (
const size_t ic) {
420 ordinal_type eOrt[12];
421 orts(ic).getEdgeOrientation(eOrt, numEdges);
422 ordinal_type ort = eOrt[ie];
425 edgeEPoints, subcellParamEdge, topoKey, ie, ort);
429 for(ordinal_type iface=0; iface<numFaces; ++iface) {
430 auto facePointsRange = ePointsRange(faceDim, iface);
431 auto faceEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getEvalPoints(faceDim,iface,ePointType));
433 const auto topoKey = refTopologyKey(faceDim,iface);
436 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
437 KOKKOS_LAMBDA (
const size_t ic) {
438 ordinal_type fOrt[6];
439 orts(ic).getFaceOrientation(fOrt, numFaces);
440 ordinal_type ort = fOrt[iface];
443 faceEPoints, subcellParamFace, topoKey, iface, ort);
449 auto pointsRange = ePointsRange(dim, 0);
450 auto cellEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getEvalPoints(dim,0,ePointType));
454 const auto topoKey = refTopologyKey(dim,0);
458 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
459 KOKKOS_LAMBDA (
const size_t ic) {
460 ordinal_type ort = 0;
462 orts(ic).getEdgeOrientation(&ort,1);
464 orts(ic).getFaceOrientation(&ort,1);
478 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtTargetEPoints,
479 const typename BasisType::ScalarViewType targetEPoints,
480 const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
481 const BasisType* cellBasis,
484 typedef typename BasisType::scalarType
scalarType;
485 typedef Kokkos::DynRankView<scalarType,DeviceType> ScalarViewType;
486 typedef Kokkos::pair<ordinal_type,ordinal_type> range_type;
488 ordinal_type dim = cellTopo.getDimension();
489 ordinal_type numTotalTargetEPoints(targetAtTargetEPoints.extent(1));
490 ordinal_type basisCardinality = cellBasis->getCardinality();
491 ordinal_type numCells = targetAtTargetEPoints.extent(0);
492 const ordinal_type edgeDim = 1;
493 const ordinal_type faceDim = 2;
494 const ordinal_type fieldDim = (targetAtTargetEPoints.rank()==2) ? 1 : targetAtTargetEPoints.extent(2);
496 ordinal_type numVertices = (cellBasis->getDofCount(0, 0) > 0) ? cellTopo.getVertexCount() : 0;
497 ordinal_type numEdges = (cellBasis->getDofCount(1, 0) > 0) ? cellTopo.getEdgeCount() : 0;
498 ordinal_type numFaces = (cellBasis->getDofCount(2, 0) > 0) ? cellTopo.getFaceCount() : 0;
500 ScalarViewType refEdgesVec(
"refEdgesVec", numEdges, dim);
501 ScalarViewType refFacesTangents(
"refFaceTangents", numFaces, dim, 2);
502 ScalarViewType refFacesNormal(
"refFaceNormal", numFaces, dim);
504 ordinal_type numVertexDofs = numVertices;
506 ordinal_type numEdgeDofs(0);
507 for(ordinal_type ie=0; ie<numEdges; ++ie)
508 numEdgeDofs += cellBasis->getDofCount(edgeDim,ie);
510 ordinal_type numFaceDofs(0);
511 for(ordinal_type iface=0; iface<numFaces; ++iface)
512 numFaceDofs += cellBasis->getDofCount(faceDim,iface);
514 Kokkos::View<ordinal_type*, DeviceType> computedDofs(
"computedDofs", numVertexDofs+numEdgeDofs+numFaceDofs);
521 ScalarViewType basisEPoints(
"basisEPoints",numCells,numTotalBasisEPoints, dim);
524 auto tagToOrdinal = Kokkos::create_mirror_view_and_copy(MemSpaceType(), cellBasis->getAllDofOrdinal());
526 ScalarViewType basisAtBasisEPoints(
"basisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints, fieldDim);
527 ScalarViewType basisAtTargetEPoints(
"basisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints, fieldDim);
530 ScalarViewType nonOrientedBasisAtBasisEPoints(
"nonOrientedBasisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints);
531 ScalarViewType nonOrientedBasisAtTargetEPoints(
"nonOrientedBasisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints);
532 for(ordinal_type ic=0; ic<numCells; ++ic) {
533 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtTargetEPoints,ic,Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(targetEPoints, ic, Kokkos::ALL(), Kokkos::ALL()));
534 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtBasisEPoints,ic,Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(basisEPoints, ic, Kokkos::ALL(), Kokkos::ALL()));
537 Kokkos::ALL(),0), nonOrientedBasisAtBasisEPoints, orts, cellBasis);
539 Kokkos::ALL(), Kokkos::ALL(),0), nonOrientedBasisAtTargetEPoints, orts, cellBasis);
542 ScalarViewType nonOrientedBasisAtBasisEPoints(
"nonOrientedBasisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints,fieldDim);
543 ScalarViewType nonOrientedBasisAtTargetEPoints(
"nonOrientedBasisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints,fieldDim);
544 for(ordinal_type ic=0; ic<numCells; ++ic) {
545 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtTargetEPoints,ic,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(targetEPoints, ic, Kokkos::ALL(), Kokkos::ALL()));
546 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtBasisEPoints,ic,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), Kokkos::subview(basisEPoints, ic, Kokkos::ALL(), Kokkos::ALL()));
554 auto hostComputedDofs = Kokkos::create_mirror_view(computedDofs);
555 ordinal_type computedDofsCount = 0;
556 for(ordinal_type iv=0; iv<numVertices; ++iv)
557 hostComputedDofs(computedDofsCount++) = cellBasis->getDofOrdinal(0, iv, 0);
559 for(ordinal_type ie=0; ie<numEdges; ++ie) {
560 ordinal_type edgeCardinality = cellBasis->getDofCount(edgeDim,ie);
561 for(ordinal_type i=0; i<edgeCardinality; ++i)
562 hostComputedDofs(computedDofsCount++) = cellBasis->getDofOrdinal(edgeDim, ie, i);
565 for(ordinal_type iface=0; iface<numFaces; ++iface) {
566 ordinal_type faceCardinality = cellBasis->getDofCount(faceDim,iface);
567 for(ordinal_type i=0; i<faceCardinality; ++i)
568 hostComputedDofs(computedDofsCount++) = cellBasis->getDofOrdinal(faceDim, iface, i);
570 Kokkos::deep_copy(computedDofs,hostComputedDofs);
573 bool isHGradBasis = (cellBasis->getFunctionSpace() == FUNCTION_SPACE_HGRAD);
574 bool isHCurlBasis = (cellBasis->getFunctionSpace() == FUNCTION_SPACE_HCURL);
575 bool isHDivBasis = (cellBasis->getFunctionSpace() == FUNCTION_SPACE_HDIV);
576 ordinal_type faceDofDim = isHCurlBasis ? 3 : 1;
577 ScalarViewType edgeCoeff(
"edgeCoeff", fieldDim);
580 const Kokkos::RangePolicy<ExecSpaceType> policy(0, numCells);
584 auto targetEPointsRange = Kokkos::create_mirror_view_and_copy(MemSpaceType(), projStruct->
getTargetPointsRange());
586 decltype(targetAtTargetEPoints),
decltype(basisAtTargetEPoints)> functorType;
587 Kokkos::parallel_for(policy, functorType(basisCoeffs, tagToOrdinal, targetEPointsRange,
588 targetAtTargetEPoints, basisAtTargetEPoints, numVertices));
592 for(ordinal_type ie=0; ie<numEdges; ++ie) {
593 auto edgeVec = Kokkos::subview(refEdgesVec, ie, Kokkos::ALL());
598 }
else if(isHDivBasis) {
601 deep_copy(edgeVec, 1.0);
604 ordinal_type edgeCardinality = cellBasis->getDofCount(edgeDim,ie);
605 ordinal_type numBasisEPoints = range_size(basisEPointsRange(edgeDim, ie));
606 ordinal_type numTargetEPoints = range_size(targetEPointsRange(edgeDim, ie));
608 ScalarViewType basisDofAtBasisEPoints(
"BasisDofAtBasisEPoints",numCells,edgeCardinality, numBasisEPoints);
609 ScalarViewType tragetDofAtTargetEPoints(
"TargetDofAtTargetEPoints",numCells, numTargetEPoints);
610 ScalarViewType weightedBasisAtBasisEPoints(
"weightedTanBasisAtBasisEPoints",numCells,edgeCardinality, numBasisEPoints);
611 ScalarViewType weightedBasisAtTargetEPoints(
"weightedTanBasisAtTargetEPoints",numCells,edgeCardinality, numTargetEPoints);
612 ScalarViewType negPartialProj(
"negPartialProj", numCells, numBasisEPoints);
614 auto targetEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getTargetEvalWeights(edgeDim,ie));
615 auto basisEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getBasisEvalWeights(edgeDim,ie));
618 ordinal_type offsetBasis = basisEPointsRange(edgeDim, ie).first;
619 ordinal_type offsetTarget = targetEPointsRange(edgeDim, ie).first;
623 decltype(computedDofs),
decltype(tagToOrdinal),
decltype(targetAtTargetEPoints)> functorTypeEdge;
625 Kokkos::parallel_for(policy, functorTypeEdge(basisCoeffs, negPartialProj, basisDofAtBasisEPoints,
626 basisAtBasisEPoints, basisEWeights, weightedBasisAtBasisEPoints, targetEWeights,
627 basisAtTargetEPoints, weightedBasisAtTargetEPoints, computedDofs, tagToOrdinal,
628 targetAtTargetEPoints,tragetDofAtTargetEPoints, refEdgesVec, fieldDim,
629 edgeCardinality, offsetBasis, offsetTarget, numVertexDofs, edgeDim, ie));
632 ScalarViewType edgeMassMat_(
"edgeMassMat_", numCells, edgeCardinality, edgeCardinality),
633 edgeRhsMat_(
"rhsMat_", numCells, edgeCardinality);
640 typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, DeviceType> WorkArrayViewType;
641 ScalarViewType t_(
"t",numCells, edgeCardinality);
642 WorkArrayViewType w_(
"w",numCells, edgeCardinality);
644 auto edgeDof = Kokkos::subview(tagToOrdinal, edgeDim, ie, Kokkos::ALL());
647 edgeSystem.
solve(basisCoeffs, edgeMassMat_, edgeRhsMat_, t_, w_, edgeDof, edgeCardinality);
650 typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParamFace;
654 for(ordinal_type iface=0; iface<numFaces; ++iface) {
655 const auto topoKey = refTopologyKey(faceDim,iface);
656 ordinal_type faceCardinality = cellBasis->getDofCount(faceDim,iface);
658 ordinal_type numTargetEPoints = range_size(targetEPointsRange(faceDim, iface));
659 ordinal_type numBasisEPoints = range_size(basisEPointsRange(faceDim, iface));
661 ScalarViewType faceBasisDofAtBasisEPoints(
"faceBasisDofAtBasisEPoints",numCells,faceCardinality, numBasisEPoints,faceDofDim);
662 ScalarViewType wBasisDofAtBasisEPoints(
"weightedBasisDofAtBasisEPoints",numCells,faceCardinality, numBasisEPoints,faceDofDim);
664 ScalarViewType faceBasisAtTargetEPoints(
"faceBasisDofAtTargetEPoints",numCells,faceCardinality, numTargetEPoints,faceDofDim);
665 ScalarViewType wBasisDofAtTargetEPoints(
"weightedBasisDofAtTargetEPoints",numCells,faceCardinality, numTargetEPoints,faceDofDim);
667 ScalarViewType targetDofAtTargetEPoints(
"targetDofAtTargetEPoints",numCells, numTargetEPoints,faceDofDim);
668 ScalarViewType negPartialProj(
"negComputedProjection", numCells,numBasisEPoints,faceDofDim);
670 ordinal_type offsetBasis = basisEPointsRange(faceDim, iface).first;
671 ordinal_type offsetTarget = targetEPointsRange(faceDim, iface).first;
672 auto targetEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getTargetEvalWeights(faceDim,iface));
673 auto basisEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getBasisEvalWeights(faceDim,iface));
677 decltype(computedDofs),
decltype(tagToOrdinal),
decltype(orts),
decltype(targetAtTargetEPoints),
decltype(subcellParamFace)> functorTypeFace;
679 Kokkos::parallel_for(policy, functorTypeFace(basisCoeffs, negPartialProj,faceBasisDofAtBasisEPoints,
680 basisAtBasisEPoints, basisEWeights, wBasisDofAtBasisEPoints, targetEWeights,
681 basisAtTargetEPoints, wBasisDofAtTargetEPoints, computedDofs, tagToOrdinal,
682 orts, targetAtTargetEPoints,targetDofAtTargetEPoints,
683 subcellParamFace, fieldDim, faceCardinality, offsetBasis,
684 offsetTarget, numVertexDofs+numEdgeDofs, numFaces, faceDim,
685 dim, iface, topoKey, isHCurlBasis, isHDivBasis));
687 typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, DeviceType> WorkArrayViewType;
688 ScalarViewType faceMassMat_(
"faceMassMat_", numCells, faceCardinality, faceCardinality),
689 faceRhsMat_(
"rhsMat_", numCells, faceCardinality);
695 ScalarViewType t_(
"t",numCells, faceCardinality);
696 WorkArrayViewType w_(
"w",numCells,faceCardinality);
698 auto faceDof = Kokkos::subview(tagToOrdinal, faceDim, iface, Kokkos::ALL());
701 faceSystem.
solve(basisCoeffs, faceMassMat_, faceRhsMat_, t_, w_, faceDof, faceCardinality);
704 ordinal_type numElemDofs = cellBasis->getDofCount(dim,0);
709 auto cellDofs = Kokkos::subview(tagToOrdinal, dim, 0, Kokkos::ALL());
711 range_type cellPointsRange = targetEPointsRange(dim, 0);
713 ordinal_type numTargetEPoints = range_size(targetEPointsRange(dim,0));
714 ordinal_type numBasisEPoints = range_size(basisEPointsRange(dim,0));
716 ScalarViewType internalBasisAtBasisEPoints(
"internalBasisAtBasisEPoints",numCells,numElemDofs, numBasisEPoints, fieldDim);
717 ScalarViewType negPartialProj(
"negPartialProj", numCells, numBasisEPoints, fieldDim);
719 auto targetEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getTargetEvalWeights(dim,0));
720 auto basisEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getBasisEvalWeights(dim,0));
721 ordinal_type offsetBasis = basisEPointsRange(dim, 0).first;
722 ordinal_type offsetTarget = targetEPointsRange(dim, 0).first;
725 ScalarViewType wBasisAtBasisEPoints(
"weightedBasisAtBasisEPoints",numCells,numElemDofs, numBasisEPoints,fieldDim);
726 ScalarViewType wBasisDofAtTargetEPoints(
"weightedBasisAtTargetEPoints",numCells,numElemDofs, numTargetEPoints,fieldDim);
728 typedef ComputeBasisCoeffsOnCells_L2<
decltype(basisCoeffs), ScalarViewType,
decltype(basisEWeights),
decltype(computedDofs),
decltype(cellDofs)> functorType;
729 Kokkos::parallel_for(policy, functorType( basisCoeffs, negPartialProj, internalBasisAtBasisEPoints,
730 basisAtBasisEPoints, basisEWeights, wBasisAtBasisEPoints, targetEWeights, basisAtTargetEPoints, wBasisDofAtTargetEPoints,
731 computedDofs, cellDofs, fieldDim, numElemDofs, offsetBasis, offsetTarget, numVertexDofs+numEdgeDofs+numFaceDofs));
733 typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, DeviceType> WorkArrayViewType;
734 ScalarViewType cellMassMat_(
"cellMassMat_", numCells, numElemDofs, numElemDofs),
735 cellRhsMat_(
"rhsMat_", numCells, numElemDofs);
740 Kokkos::subview(wBasisDofAtTargetEPoints,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),0));
745 ScalarViewType t_(
"t",numCells, numElemDofs);
746 WorkArrayViewType w_(
"w",numCells,numElemDofs);
748 cellSystem.
solve(basisCoeffs, cellMassMat_, cellRhsMat_, t_, w_, cellDofs, numElemDofs);
807 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtTargetEPoints,
808 const Kokkos::DynRankView<ortValueType, ortProperties...> orts,
809 const BasisType* cellBasis,
812 typedef typename BasisType::scalarType
scalarType;
813 typedef Kokkos::DynRankView<scalarType,DeviceType> ScalarViewType;
816 ordinal_type dim = cellTopo.getDimension();
818 auto basisEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),
820 auto targetEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),
824 ordinal_type numTotalTargetEPoints(targetAtTargetEPoints.extent(1));
825 ordinal_type basisCardinality = cellBasis->getCardinality();
826 ordinal_type numCells = targetAtTargetEPoints.extent(0);
827 const ordinal_type fieldDim = (targetAtTargetEPoints.rank()==2) ? 1 : targetAtTargetEPoints.extent(2);
830 ScalarViewType basisAtBasisEPoints(
"basisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints, fieldDim);
831 ScalarViewType basisAtTargetEPoints(
"basisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints, fieldDim);
834 ScalarViewType nonOrientedBasisAtBasisEPoints(
"nonOrientedBasisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints);
835 ScalarViewType nonOrientedBasisAtTargetEPoints(
"nonOrientedBasisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints);
836 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtTargetEPoints,0,Kokkos::ALL(),Kokkos::ALL()), targetEPoints);
837 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtBasisEPoints,0,Kokkos::ALL(),Kokkos::ALL()), basisEPoints);
842 Kokkos::ALL(),0), nonOrientedBasisAtBasisEPoints, orts, cellBasis);
844 Kokkos::ALL(), Kokkos::ALL(),0), nonOrientedBasisAtTargetEPoints, orts, cellBasis);
847 ScalarViewType nonOrientedBasisAtBasisEPoints(
"nonOrientedBasisAtBasisEPoints",numCells,basisCardinality, numTotalBasisEPoints,fieldDim);
848 ScalarViewType nonOrientedBasisAtTargetEPoints(
"nonOrientedBasisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints,fieldDim);
849 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtTargetEPoints,0,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), targetEPoints);
850 cellBasis->getValues(Kokkos::subview(nonOrientedBasisAtBasisEPoints,0,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), basisEPoints);
852 RealSpaceTools<DeviceType>::clone(nonOrientedBasisAtTargetEPoints, Kokkos::subview(nonOrientedBasisAtTargetEPoints,0,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()));
859 const Kokkos::RangePolicy<ExecSpaceType> policy(0, numCells);
860 ordinal_type numElemDofs = cellBasis->getCardinality();
862 auto targetEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getTargetEvalWeights(dim,0));
863 auto basisEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getBasisEvalWeights(dim,0));
865 ScalarViewType wBasisAtBasisEPoints(
"weightedBasisAtBasisEPoints",numCells,numElemDofs, numTotalBasisEPoints,fieldDim);
866 ScalarViewType wBasisDofAtTargetEPoints(
"weightedBasisAtTargetEPoints",numCells,numElemDofs, numTotalTargetEPoints,fieldDim);
869 Kokkos::parallel_for(
"Multiply basis by weights", policy, functorType(basisAtBasisEPoints, basisEWeights,
870 wBasisAtBasisEPoints, targetEWeights, basisAtTargetEPoints, wBasisDofAtTargetEPoints, fieldDim, numElemDofs));
872 typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, DeviceType> WorkArrayViewType;
873 ScalarViewType cellMassMat_(
"cellMassMat_", numCells, numElemDofs, numElemDofs),
874 cellRhsMat_(
"rhsMat_", numCells, numElemDofs);
879 Kokkos::subview(wBasisDofAtTargetEPoints,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),0));
883 Kokkos::DynRankView<ordinal_type,Kokkos::HostSpace> hCellDofs(
"cellDoFs", numElemDofs);
884 for(ordinal_type i=0; i<numElemDofs; ++i) hCellDofs(i) = i;
885 auto cellDofs = Kokkos::create_mirror_view_and_copy(MemSpaceType(),hCellDofs);
887 ScalarViewType t_(
"t",numCells, numElemDofs);
888 WorkArrayViewType w_(
"w",numCells,numElemDofs);
890 cellSystem.
solve(basisCoeffs, cellMassMat_, cellRhsMat_, t_, w_, cellDofs, numElemDofs);
897 const targetViewType targetAtTargetEPoints,
898 const BasisType* cellBasis,
901 typedef typename BasisType::scalarType
scalarType;
902 typedef Kokkos::DynRankView<scalarType,DeviceType> ScalarViewType;
905 ordinal_type dim = cellTopo.getDimension();
907 auto basisEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),
909 auto targetEPoints = Kokkos::create_mirror_view_and_copy(MemSpaceType(),
912 ordinal_type numTotalTargetEPoints(targetAtTargetEPoints.extent(1));
913 ordinal_type basisCardinality = cellBasis->getCardinality();
914 ordinal_type numCells = targetAtTargetEPoints.extent(0);
915 const ordinal_type fieldDim = (targetAtTargetEPoints.rank()==2) ? 1 : targetAtTargetEPoints.extent(2);
918 ScalarViewType basisAtBasisEPoints(
"basisAtBasisEPoints",1,basisCardinality, numTotalBasisEPoints, fieldDim);
919 ScalarViewType basisAtTargetEPoints(
"basisAtTargetEPoints",numCells,basisCardinality, numTotalTargetEPoints, fieldDim);
922 cellBasis->getValues(Kokkos::subview(basisAtTargetEPoints,0,Kokkos::ALL(),Kokkos::ALL(),0), targetEPoints);
923 cellBasis->getValues(Kokkos::subview(basisAtBasisEPoints,0,Kokkos::ALL(),Kokkos::ALL(),0), basisEPoints);
928 cellBasis->getValues(Kokkos::subview(basisAtTargetEPoints,0,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), targetEPoints);
929 cellBasis->getValues(Kokkos::subview(basisAtBasisEPoints,0,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL()), basisEPoints);
935 const Kokkos::RangePolicy<ExecSpaceType> policy(0, numCells);
936 ordinal_type numElemDofs = cellBasis->getCardinality();
938 auto targetEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getTargetEvalWeights(dim,0));
939 auto basisEWeights = Kokkos::create_mirror_view_and_copy(MemSpaceType(),projStruct->
getBasisEvalWeights(dim,0));
941 ScalarViewType wBasisAtBasisEPoints(
"weightedBasisAtBasisEPoints",1,numElemDofs, numTotalBasisEPoints,fieldDim);
942 ScalarViewType wBasisDofAtTargetEPoints(
"weightedBasisAtTargetEPoints",numCells,numElemDofs, numTotalTargetEPoints,fieldDim);
943 Kokkos::DynRankView<ordinal_type, DeviceType> cellDofs(
"cellDoFs", numElemDofs);
945 Kokkos::parallel_for(Kokkos::RangePolicy<ExecSpaceType>(0,numElemDofs),
946 KOKKOS_LAMBDA (
const int &j) {
947 for(ordinal_type d=0; d <fieldDim; ++d) {
948 for(ordinal_type iq=0; iq < numTotalBasisEPoints; ++iq)
949 wBasisAtBasisEPoints(0,j,iq,d) = basisAtBasisEPoints(0,j,iq,d) * basisEWeights(iq);
950 for(ordinal_type iq=0; iq <numTotalTargetEPoints; ++iq) {
951 wBasisDofAtTargetEPoints(0,j,iq,d) = basisAtTargetEPoints(0,j,iq,d)* targetEWeights(iq);
958 typedef Kokkos::DynRankView<scalarType, Kokkos::LayoutRight, DeviceType> WorkArrayViewType;
959 ScalarViewType cellMassMat_(
"cellMassMat_", 1, numElemDofs, numElemDofs),
960 cellRhsMat_(
"rhsMat_", numCells, numElemDofs);
965 Kokkos::subview(wBasisDofAtTargetEPoints,Kokkos::ALL(),Kokkos::ALL(),Kokkos::ALL(),0));
969 ScalarViewType t_(
"t",1, numElemDofs);
970 WorkArrayViewType w_(
"w",numCells,numElemDofs);
972 cellSystem.
solve(basisCoeffs, cellMassMat_, cellRhsMat_, t_, w_, cellDofs, numElemDofs);