Intrepid2
Intrepid2_ProjectionTools.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
14#ifndef __INTREPID2_PROJECTIONTOOLS_HPP__
15#define __INTREPID2_PROJECTIONTOOLS_HPP__
16
17#include "Intrepid2_ConfigDefs.hpp"
18#include "Intrepid2_Types.hpp"
19#include "Intrepid2_Utils.hpp"
20
21#include "Shards_CellTopology.hpp"
22#include "Shards_BasicTopologies.hpp"
23
25
26#include "Intrepid2_Basis.hpp"
27
29
30// -- Lower order family
33
36
40
44
45#include "Teuchos_LAPACK.hpp"
47
49
50#ifdef HAVE_INTREPID2_KOKKOSKERNELS
51#include "KokkosBatched_QR_Serial_Internal.hpp"
52#include "KokkosBatched_ApplyQ_Serial_Internal.hpp"
53#include "KokkosBatched_Trsv_Decl.hpp"
54#include "KokkosBatched_Util.hpp"
55#endif
56
57namespace Intrepid2 {
58
111template<typename DeviceType>
113public:
114 using ExecSpaceType = typename DeviceType::execution_space;
115 using MemSpaceType = typename DeviceType::memory_space;
116 using EvalPointsType = typename ProjectionStruct<DeviceType, double>::EvalPointsType;
117
118
138 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
139 typename funValsValueType, class ...funValsProperties,
140 typename BasisType,
141 typename ortValueType, class ...ortProperties>
142 static void
143 getL2BasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
144 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
145 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
146 const BasisType* cellBasis,
148
149
173 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
174 typename funValsValueType, class ...funValsProperties,
175 typename BasisType,
176 typename ortValueType, class ...ortProperties>
177 static void
178 getL2DGBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
179 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
180 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
181 const BasisType* cellBasis,
183
184
207 template<typename basisViewType, typename targetViewType, typename BasisType>
208 static void
209 getL2DGBasisCoeffs(basisViewType basisCoeffs,
210 const targetViewType targetAtTargetEPoints,
211 const BasisType* cellBasis,
213
235 template<class BasisCoeffsViewType, class TargetValueViewType, class TargetGradViewType,
236 class BasisType, class OrientationViewType>
237 static void
238 getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs,
239 const TargetValueViewType targetAtEvalPoints,
240 const TargetGradViewType targetGradAtGradEvalPoints,
241 const OrientationViewType cellOrientations,
242 const BasisType* cellBasis,
244
245
269 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
270 typename funValsValueType, class ...funValsProperties,
271 typename BasisType,
272 typename ortValueType, class ...ortProperties>
273 static void
274 getHCurlBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
275 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
276 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetCurlAtCurlEvalPoints,
277 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
278 const BasisType* cellBasis,
280
302 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
303 typename funValsValueType, class ...funValsProperties,
304 typename BasisType,
305 typename ortValueType, class ...ortProperties>
306 static void
307 getHDivBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
308 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
309 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetDivAtDivEvalPoints,
310 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
311 const BasisType* cellBasis,
313
314
332 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
333 typename funValsValueType, class ...funValsProperties,
334 typename BasisType,
335 typename ortValueType, class ...ortProperties>
336 static void
337 getHVolBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
338 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
339 [[maybe_unused]] const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
340 const BasisType* cellBasis,
342
343
344
360 template<typename dstViewType,
361 typename dstBasisType,
362 typename srcViewType,
363 typename srcBasisType,
364 typename OrientationViewType>
365 static void
366 projectField(dstViewType dstCoeffs,
367 const dstBasisType* dstBasis,
368 const srcViewType srcCoeffs,
369 const srcBasisType* srcBasis,
370 const OrientationViewType cellOrientations){
371
372
373 INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getFunctionSpace() != srcBasis->getFunctionSpace(), std::runtime_error,
374 "The source and destination bases are not compatible. They need to belong to the same function space");
375 INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getBaseCellTopology().getKey() != srcBasis->getBaseCellTopology().getKey(), std::runtime_error,
376 "The source and destination bases are not compatible. They do not have the same basic cell topology");
377
379 projStruct.createL2ProjectionStruct(dstBasis, srcBasis->getDegree());
380
381
382 ordinal_type numCells = cellOrientations.extent(0);
383 ordinal_type dim = srcBasis->getBaseCellTopology().getDimension();
384 ordinal_type srcBasisCardinality = srcBasis->getCardinality();
385 ordinal_type fieldDimension = (srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HCURL || srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HDIV) ? dim : 1;
386
387 auto evaluationPoints = projStruct.getAllEvalPoints();
388 ordinal_type numPoints = evaluationPoints.extent(0);
389
390 using outViewType = Kokkos::DynRankView<typename srcBasisType::OutputValueType, DeviceType>;
391 outViewType srcAtEvalPoints, refBasisAtEvalPoints, basisAtEvalPoints;
392 if(fieldDimension == dim) {
393 srcAtEvalPoints = outViewType("srcAtEvalPoints", numCells, numPoints, dim);
394 refBasisAtEvalPoints = outViewType("refBasisAtEvalPoints", srcBasisCardinality, numPoints, dim);
395 basisAtEvalPoints = outViewType("basisAtEvalPoints", numCells, srcBasisCardinality, numPoints, dim);
396 } else {
397 srcAtEvalPoints = outViewType("srcAtEvalPoints", numCells, numPoints);
398 refBasisAtEvalPoints = outViewType("refBasisAtEvalPoints", srcBasisCardinality, numPoints);
399 basisAtEvalPoints = outViewType("basisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
400 }
401
402 srcBasis->getValues(refBasisAtEvalPoints,evaluationPoints);
403
404 // Modify basis values to account for orientations
406 refBasisAtEvalPoints,
407 cellOrientations,
408 srcBasis);
409
410 Kokkos::parallel_for(Kokkos::RangePolicy<typename DeviceType::execution_space>(0,numCells),
411 KOKKOS_LAMBDA (const int &ic) {
412 for(int j=0; j<numPoints; ++j) {
413 for(int k=0; k<srcBasisCardinality; ++k) {
414 for(int d=0; d<fieldDimension; ++d)
415 srcAtEvalPoints.access(ic,j,d) += srcCoeffs(ic,k)*basisAtEvalPoints.access(ic,k,j,d);
416 }
417 }
418 });
419 ExecSpaceType().fence();
420
421 getL2BasisCoeffs(dstCoeffs,
422 srcAtEvalPoints,
423 cellOrientations,
424 dstBasis,
425 &projStruct);
426 }
427
428
429
439 struct ElemSystem {
440
441
442 std::string systemName_;
443 bool matrixIndependentOfCell_;
444
452 ElemSystem (std::string systemName, bool matrixIndependentOfCell) :
453 systemName_(systemName), matrixIndependentOfCell_(matrixIndependentOfCell){};
454
455
456
482 template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
483 void solve(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 tau,
484 ViewType3 w,const ViewType4 elemDof, ordinal_type n, ordinal_type m=0) {
485#ifdef HAVE_INTREPID2_KOKKOSKERNELS
486 solveDevice(basisCoeffs, elemMat, elemRhs, tau,
487 w, elemDof, n, m);
488#else
489 solveHost(basisCoeffs, elemMat, elemRhs, tau,
490 w, elemDof, n, m);
491#endif
492
493 }
494
497#ifdef HAVE_INTREPID2_KOKKOSKERNELS
498 template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
499 void solveDevice(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 taul,
500 ViewType3 work,const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
501 using HostDeviceType = Kokkos::Device<Kokkos::DefaultHostExecutionSpace,Kokkos::HostSpace>;
502
503 ordinal_type numCells = basisCoeffs.extent(0);
504
505 if(matrixIndependentOfCell_) {
506 auto A0 = Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL());
507 auto tau0 = Kokkos::subview(taul, 0, Kokkos::ALL());
508
509 Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> A0_host("A0_host", A0.extent(0),A0.extent(1));
510 auto A0_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), A0_host);
511 Kokkos::deep_copy(A0_device, A0);
512 Kokkos::deep_copy(A0_host, A0_device);
513
514 for(ordinal_type i=n; i<n+m; ++i)
515 for(ordinal_type j=0; j<n; ++j)
516 A0_host(i,j) = A0_host(j,i);
517
518 Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> tau0_host("A0_host", tau0.extent(0));
519 auto tau0_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), tau0_host);
520 auto w0_host = Kokkos::create_mirror_view(Kokkos::subview(work, 0, Kokkos::ALL()));
521
522 //computing QR of A0. QR is saved in A0 and tau0
523 KokkosBatched::Impl::SerialQR_Internal::invoke
524 (A0_host.extent(0), A0_host.extent(1),
525 A0_host.data(), A0_host.stride(0), A0_host.stride(1),
526 tau0_host.data(), tau0_host.stride(0), w0_host.data());
527
528 Kokkos::deep_copy(A0_device, A0_host);
529 Kokkos::deep_copy(A0, A0_device);
530 Kokkos::deep_copy(tau0_device, tau0_host);
531 Kokkos::deep_copy(tau0, tau0_device);
532
533 Kokkos::parallel_for (systemName_,
534 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
535 KOKKOS_LAMBDA (const size_t ic) {
536 auto w = Kokkos::subview(work, ic, Kokkos::ALL());
537
538 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
539
540 //b'*Q0 -> b
541 KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
542 1, A0.extent(0), A0.extent(1),
543 A0.data(), A0.stride(0), A0.stride(1),
544 tau0.data(), tau0.stride(0),
545 b.data(), 1, b.stride(0),
546 w.data());
547
548 // R0^{-1} b -> b
549 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(1.0, A0, b);
550
551 //scattering b into the basis coefficients
552 for(ordinal_type i=0; i<n; ++i){
553 basisCoeffs(ic,elemDof(i)) = b(i);
554 }
555 });
556
557 } else {
558
559 Kokkos::parallel_for (systemName_,
560 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
561 KOKKOS_LAMBDA (const size_t ic) {
562
563 auto A = Kokkos::subview(elemMat, ic, Kokkos::ALL(), Kokkos::ALL());
564 auto tau = Kokkos::subview(taul, ic, Kokkos::ALL());
565 auto w = Kokkos::subview(work, ic, Kokkos::ALL());
566
567 for(ordinal_type i=n; i<n+m; ++i)
568 for(ordinal_type j=0; j<n; ++j)
569 A(i,j) = A(j,i);
570
571 //computing QR of A. QR is saved in A and tau
572 KokkosBatched::Impl::SerialQR_Internal::invoke
573 (A.extent(0), A.extent(1),
574 A.data(), A.stride(0), A.stride(1), tau.data(), tau.stride(0), w.data());
575
576 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
577
578 //b'*Q -> b
579 KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
580 1, A.extent(0), A.extent(1),
581 A.data(), A.stride(0), A.stride(1),
582 tau.data(), tau.stride(0),
583 b.data(), 1, b.stride(0),
584 w.data());
585
586 // R^{-1} b -> b
587 KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(1.0, A, b);
588
589 //scattering b into the basis coefficients
590 for(ordinal_type i=0; i<n; ++i){
591 basisCoeffs(ic,elemDof(i)) = b(i);
592 }
593 });
594 }
595 }
596#endif
597
601 template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
602 void solveHost(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 ,
603 ViewType3, const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
604 using value_type = typename ViewType2::value_type;
605 using device_type = DeviceType;
606 using host_exec_space = Kokkos::DefaultHostExecutionSpace;
607 using host_memory_space = Kokkos::HostSpace;
608 using host_device_type = Kokkos::Device<host_exec_space,host_memory_space>;
609 using vector_host_type = Kokkos::View<value_type*,host_device_type>;
610 using scratch_host_type = Kokkos::View<value_type*,host_exec_space::scratch_memory_space>;
611 using matrix_host_type = Kokkos::View<value_type**,Kokkos::LayoutLeft,host_device_type>;
612 using do_not_init_tag = Kokkos::ViewAllocateWithoutInitializing;
613 using host_team_policy_type = Kokkos::TeamPolicy<host_exec_space>;
614 using host_range_policy_type = Kokkos::RangePolicy<host_exec_space>;
615
617 Kokkos::fence();
618
620 const ordinal_type numCells = basisCoeffs.extent(0);
621 const ordinal_type numRows = m+n, numCols = n;
622
624 Teuchos::LAPACK<ordinal_type,value_type> lapack;
625
627 Kokkos::View<ordinal_type*,host_device_type> elemDof_host(do_not_init_tag("elemDof_host"), elemDof.extent(0));
628 {
629 auto elemDof_device = Kokkos::create_mirror_view(typename device_type::memory_space(), elemDof_host);
630 Kokkos::deep_copy(elemDof_device, elemDof); Kokkos::fence();
631 Kokkos::deep_copy(elemDof_host, elemDof_device);
632 }
633
635 auto elemRhs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemRhs);
636 auto elemMat_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemMat);
637
639 auto basisCoeffs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), basisCoeffs);
640
641 if (matrixIndependentOfCell_) {
643 matrix_host_type A(do_not_init_tag("A"), numRows, numRows);
644 {
645 for (ordinal_type j=0;j<numRows;++j)
646 for (ordinal_type i=0;i<numRows;++i)
647 A(i, j) = elemMat_host(0, i, j);
648
649 for (ordinal_type j=0;j<numCols;++j)
650 for (ordinal_type i=numCols;i<numRows;++i)
651 A(i, j) = A(j, i);
652 }
653
654 ordinal_type lwork(-1);
655 {
656 ordinal_type info(0);
657 value_type work[2];
658 lapack.GELS('N',
659 numRows, numRows, numCells,
660 nullptr, std::max(1,numRows),
661 nullptr, std::max(1,numRows),
662 &work[0], lwork,
663 &info);
664 lwork = work[0];
665 }
666
667 matrix_host_type C(do_not_init_tag("C"), numRows, numCells);
668
669 host_range_policy_type policy(0, numCells);
670 {
671 Kokkos::parallel_for
672 ("ProjectionTools::solveHost::matrixIndependentOfCell::pack",
673 policy, [=](const ordinal_type & ic) {
674 for (ordinal_type i=0;i<numRows;++i)
675 C(i,ic) = elemRhs_host(ic, i);
676 });
677 }
678 {
680 vector_host_type work(do_not_init_tag("work"), lwork);
681 ordinal_type info(0);
682 lapack.GELS('N',
683 numRows, numRows, numCells,
684 A.data(), std::max(1,numRows),
685 C.data(), std::max(1,numRows),
686 work.data(), lwork,
687 &info);
688 INTREPID2_TEST_FOR_EXCEPTION
689 (info != 0, std::runtime_error, "GELS return non-zero info code");
690 }
691 {
692 Kokkos::parallel_for
693 ("ProjectionTools::solveHost::matrixIndependentOfCell::unpack",
694 policy, [=](const ordinal_type & ic) {
695 for (ordinal_type i=0;i<numCols;++i)
696 basisCoeffs_host(ic,elemDof_host(i)) = C(i,ic);
697 });
698 }
699 } else {
700 const ordinal_type level(0);
701 host_team_policy_type policy(numCells, 1, 1);
702
704 ordinal_type lwork(-1);
705 {
706 ordinal_type info(0);
707 value_type work[2];
708 lapack.GELS('N',
709 numRows, numRows, 1,
710 nullptr, std::max(1,numRows),
711 nullptr, std::max(1,numRows),
712 &work[0], lwork,
713 &info);
714 lwork = work[0];
715 }
716
717 const ordinal_type per_team_extent = numRows*numRows + numRows + lwork;
718 const ordinal_type per_team_scratch = scratch_host_type::shmem_size(per_team_extent);
719 policy.set_scratch_size(level, Kokkos::PerTeam(per_team_scratch));
720
722 Kokkos::parallel_for
723 ("ProjectionTools::solveHost::matrixDependentOfCell",
724 policy, [=](const typename host_team_policy_type::member_type& member) {
725 const ordinal_type ic = member.league_rank();
726
727 scratch_host_type scratch(member.team_scratch(level), per_team_extent);
728 value_type * sptr = scratch.data();
729
731 matrix_host_type A(sptr, numRows, numRows); sptr += A.span();
732 for (ordinal_type j=0;j<numRows;++j)
733 for (ordinal_type i=0;i<numRows;++i)
734 A(i, j) = elemMat_host(ic, i, j);
735
736 for (ordinal_type j=0;j<numCols;++j)
737 for (ordinal_type i=numCols;i<numRows;++i)
738 A(i, j) = A(j, i);
739
740 vector_host_type c(sptr, numRows); sptr += c.span();
741 for (ordinal_type i=0;i<numRows;++i)
742 c(i) = elemRhs_host(ic, i);
743
744 vector_host_type work(sptr, lwork); sptr += work.span();
745 ordinal_type info(0);
746 lapack.GELS('N',
747 numRows, numRows, 1,
748 A.data(), std::max(1,numRows),
749 c.data(), std::max(1,numRows),
750 work.data(), lwork,
751 &info);
752 INTREPID2_TEST_FOR_EXCEPTION
753 (info != 0, std::runtime_error, "GELS return non-zero info code");
754
756 for (ordinal_type i=0;i<numCols;++i)
757 basisCoeffs_host(ic,elemDof_host(i)) = c(i);
758 });
759 }
760 Kokkos::deep_copy(basisCoeffs, basisCoeffs_host);
761 }
762 };
763
764};
765
766} //Intrepid2
767
768
769// include templated function definitions
775
776#endif
777
778
779
780
781
Header file for the abstract base class Intrepid2::Basis.
Header file for the Intrepid2::Basis_HCURL_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_WEDGE_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_WEDGE_I1_FEM class.
Stateless class that acts as a factory for a family of nodal bases (hypercube topologies only at this...
Header file for the Intrepid2::OrientationTools and Intrepid2::Impl::OrientationTools classes.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Header file for the Intrepid2::ProjectionStruct.
Header file for the Intrepid2::ProjectionTools containing definitions for HCURL projections.
Header file for the Intrepid2::ProjectionTools containing definitions for HDIV projections.
Header file for the Intrepid2::ProjectionTools containing definitions for HGRAD projections.
Header file for the Intrepid2::ProjectionTools containing definitions for HVOL projections.
Header file for the Intrepid2::ProjectionTools containing definitions for L2 projections.
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
static void modifyBasisByOrientation(Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis, const bool transpose=false)
Modify basis due to orientation.
An helper class to compute the evaluation points and weights needed for performing projections.
void createL2ProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
Initialize the ProjectionStruct for L2 projections.
view_type getAllEvalPoints(EvalPointsType type=TARGET) const
Returns the basis/target evaluation points in the cell.
A class providing static members to perform projection-based interpolations:
static void projectField(dstViewType dstCoeffs, const dstBasisType *dstBasis, const srcViewType srcCoeffs, const srcBasisType *srcBasis, const OrientationViewType cellOrientations)
Computes the L2 projection of a finite element field into a compatible finite element space.
static void getL2BasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the L2 projection of the target function.
static void getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs, const TargetValueViewType targetAtEvalPoints, const TargetGradViewType targetGradAtGradEvalPoints, const OrientationViewType cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HGrad projection of the target function.
static void getL2DGBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes evaluation points for local L2 projection for broken HGRAD HCURL HDIV and HVOL spaces.
static void getHDivBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetDivAtDivEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HDiv projection of the target function.
static void getHCurlBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetCurlAtCurlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HCurl projection of the target function.
static void getHVolBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HVol projection of the target function.
Class to solve a square system A x = b on each cell A is expected to be saddle a point (KKT) matrix o...
ElemSystem(std::string systemName, bool matrixIndependentOfCell)
Functor constructor.
void solve(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 tau, ViewType3 w, const ViewType4 elemDof, ordinal_type n, ordinal_type m=0)
Solve the system and returns the basis coefficients solve the system either using Kokkos Kernel QR or...
void solveHost(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2, ViewType3, const ViewType4 elemDof, ordinal_type n, ordinal_type m)
Parallel implementation of solve, using Kokkos Kernels QR factoriation.