47#ifndef __INTREPID2_PROJECTIONTOOLS_HPP__
48#define __INTREPID2_PROJECTIONTOOLS_HPP__
50#include "Intrepid2_ConfigDefs.hpp"
54#include "Shards_CellTopology.hpp"
55#include "Shards_BasicTopologies.hpp"
100#include "Teuchos_LAPACK.hpp"
105#ifdef HAVE_INTREPID2_KOKKOSKERNELS
106#include "KokkosBatched_QR_Serial_Internal.hpp"
107#include "KokkosBatched_ApplyQ_Serial_Internal.hpp"
108#include "KokkosBatched_Trsv_Serial_Internal.hpp"
109#include "KokkosBatched_Util.hpp"
114namespace Experimental {
182template<
typename DeviceType>
185 using ExecSpaceType =
typename DeviceType::execution_space;
186 using MemSpaceType =
typename DeviceType::memory_space;
187 using EvalPointsType =
typename ProjectionStruct<DeviceType, double>::EvalPointsType;
206 template<
typename BasisType,
207 typename ortValueType,
class ...ortProperties>
210 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
211 const BasisType* cellBasis,
213 const EvalPointsType evalPointType = EvalPointsType::TARGET
235 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
236 typename funValsValueType,
class ...funValsProperties,
238 typename ortValueType,
class ...ortProperties>
240 getL2BasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
241 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
242 const typename BasisType::ScalarViewType evaluationPoints,
243 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
244 const BasisType* cellBasis,
264 template<
typename BasisType>
267 const BasisType* cellBasis,
269 const EvalPointsType evalPointType = EvalPointsType::TARGET
295 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
296 typename funValsValueType,
class ...funValsProperties,
298 typename ortValueType,
class ...ortProperties>
300 getL2DGBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
301 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
302 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
303 const BasisType* cellBasis,
328 template<
typename basisViewType,
typename targetViewType,
typename BasisType>
331 const targetViewType targetAtTargetEPoints,
332 const BasisType* cellBasis,
355 template<
typename BasisType,
typename OrientationViewType >
358 typename BasisType::ScalarViewType gradEvalPoints,
359 const OrientationViewType cellOrientations,
360 const BasisType* cellBasis,
362 const EvalPointsType evalPointType = EvalPointsType::TARGET
388 template<
class BasisCoeffsViewType,
class TargetValueViewType,
class TargetGradViewType,
389 class BasisType,
class OrientationViewType>
392 const TargetValueViewType targetAtEvalPoints,
393 const TargetGradViewType targetGradAtGradEvalPoints,
394 const typename BasisType::ScalarViewType evaluationPoints,
395 const typename BasisType::ScalarViewType gradEvalPoints,
396 const OrientationViewType cellOrientations,
397 const BasisType* cellBasis,
420 template<
typename BasisType,
421 typename ortValueType,
class ...ortProperties>
424 typename BasisType::ScalarViewType curlEvalPoints,
425 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
426 const BasisType* cellBasis,
428 const EvalPointsType evalPointType = EvalPointsType::TARGET
456 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
457 typename funValsValueType,
class ...funValsProperties,
459 typename ortValueType,
class ...ortProperties>
461 getHCurlBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
462 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
463 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetCurlAtCurlEvalPoints,
464 const typename BasisType::ScalarViewType evaluationPoints,
465 const typename BasisType::ScalarViewType curlEvalPoints,
466 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
467 const BasisType* cellBasis,
490 template<
typename BasisType,
491 typename ortValueType,
class ...ortProperties>
494 typename BasisType::ScalarViewType divEvalPoints,
495 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
496 const BasisType* cellBasis,
498 const EvalPointsType evalPointType = EvalPointsType::TARGET
524 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
525 typename funValsValueType,
class ...funValsProperties,
527 typename ortValueType,
class ...ortProperties>
529 getHDivBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
530 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
531 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetDivAtDivEvalPoints,
532 const typename BasisType::ScalarViewType evaluationPoints,
533 const typename BasisType::ScalarViewType divEvalPoints,
534 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
535 const BasisType* cellBasis,
554 template<
typename BasisType,
555 typename ortValueType,
class ...ortProperties>
558 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
559 const BasisType* cellBasis,
561 const EvalPointsType evalPointType = EvalPointsType::TARGET
582 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
583 typename funValsValueType,
class ...funValsProperties,
585 typename ortValueType,
class ...ortProperties>
587 getHVolBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
588 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
589 const typename BasisType::ScalarViewType evaluationPoints,
590 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
591 const BasisType* cellBasis,
611 template<
typename dstViewType,
612 typename dstBasisType,
613 typename srcViewType,
614 typename srcBasisType,
615 typename OrientationViewType>
618 const dstBasisType* dstBasis,
619 const srcViewType srcCoeffs,
620 const srcBasisType* srcBasis,
621 const OrientationViewType cellOrientations){
624 INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getFunctionSpace() != srcBasis->getFunctionSpace(), std::runtime_error,
625 "The source and destination bases are not compatible. They need to belong to the same function space");
626 INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getBaseCellTopology().getKey() != srcBasis->getBaseCellTopology().getKey(), std::runtime_error,
627 "The source and destination bases are not compatible. They do not have the same basic cell topology");
633 ordinal_type numCells = cellOrientations.extent(0);
634 ordinal_type dim = srcBasis->getBaseCellTopology().getDimension();
635 ordinal_type srcBasisCardinality = srcBasis->getCardinality();
636 ordinal_type fieldDimension = (srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HCURL || srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HDIV) ? dim : 1;
638 typename Kokkos::DynRankView<typename srcBasisType::PointValueType, DeviceType> evaluationPoints(
"evaluationPoints", numCells, numPoints, dim);
644 using outViewType = Kokkos::DynRankView<typename srcBasisType::OutputValueType, DeviceType>;
645 outViewType srcAtEvalPoints, refBasisAtEvalPoints, basisAtEvalPoints;
646 if(fieldDimension == dim) {
647 srcAtEvalPoints = outViewType(
"srcAtEvalPoints", numCells, numPoints, dim);
648 refBasisAtEvalPoints = outViewType(
"refBasisAtEvalPoints", numCells, srcBasisCardinality, numPoints, dim);
649 basisAtEvalPoints = outViewType(
"basisAtEvalPoints", numCells, srcBasisCardinality, numPoints, dim);
651 srcAtEvalPoints = outViewType(
"srcAtEvalPoints", numCells, numPoints);
652 refBasisAtEvalPoints = outViewType(
"refBasisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
653 basisAtEvalPoints = outViewType(
"basisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
656 for(ordinal_type icell = 0; icell < numCells; ++icell)
657 srcBasis->getValues(Kokkos::subview(refBasisAtEvalPoints, icell, Kokkos::ALL(), Kokkos::ALL(), Kokkos::ALL()),Kokkos::subview(evaluationPoints, icell, Kokkos::ALL(), Kokkos::ALL()));
661 refBasisAtEvalPoints,
665 Kokkos::parallel_for(Kokkos::RangePolicy<typename DeviceType::execution_space>(0,numCells),
666 KOKKOS_LAMBDA (
const int &ic) {
667 for(
int j=0; j<numPoints; ++j) {
668 for(
int k=0; k<srcBasisCardinality; ++k) {
669 for(
int d=0; d<fieldDimension; ++d)
670 srcAtEvalPoints.access(ic,j,d) += srcCoeffs(ic,k)*basisAtEvalPoints.access(ic,k,j,d);
674 ExecSpaceType().fence();
698 std::string systemName_;
699 bool matrixIndependentOfCell_;
708 ElemSystem (std::string systemName,
bool matrixIndependentOfCell) :
709 systemName_(systemName), matrixIndependentOfCell_(matrixIndependentOfCell){};
738 template<
typename ViewType1,
typename ViewType2,
typename ViewType3,
typename ViewType4>
739 void solve(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 tau,
740 ViewType3 w,
const ViewType4 elemDof, ordinal_type n, ordinal_type m=0) {
741#ifdef HAVE_INTREPID2_KOKKOSKERNELS
742 solveDevice(basisCoeffs, elemMat, elemRhs, tau,
745 solveHost(basisCoeffs, elemMat, elemRhs, tau,
753#ifdef HAVE_INTREPID2_KOKKOSKERNELS
754 template<
typename ViewType1,
typename ViewType2,
typename ViewType3,
typename ViewType4>
755 void solveDevice(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 taul,
756 ViewType3 work,
const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
757 using HostDeviceType = Kokkos::Device<Kokkos::DefaultHostExecutionSpace,Kokkos::HostSpace>;
759 ordinal_type numCells = basisCoeffs.extent(0);
761 if(matrixIndependentOfCell_) {
762 auto A0 = Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL());
763 auto tau0 = Kokkos::subview(taul, 0, Kokkos::ALL());
765 Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> A0_host(
"A0_host", A0.extent(0),A0.extent(1));
766 auto A0_device = Kokkos::create_mirror_view(
typename DeviceType::memory_space(), A0_host);
767 Kokkos::deep_copy(A0_device, A0);
768 Kokkos::deep_copy(A0_host, A0_device);
770 for(ordinal_type i=n; i<n+m; ++i)
771 for(ordinal_type j=0; j<n; ++j)
772 A0_host(i,j) = A0_host(j,i);
774 Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> tau0_host(
"A0_host", tau0.extent(0));
775 auto tau0_device = Kokkos::create_mirror_view(
typename DeviceType::memory_space(), tau0_host);
776 auto w0_host = Kokkos::create_mirror_view(Kokkos::subview(work, 0, Kokkos::ALL()));
779 KokkosBatched::SerialQR_Internal::invoke(A0_host.extent(0), A0_host.extent(1),
780 A0_host.data(), A0_host.stride_0(), A0_host.stride_1(),
781 tau0_host.data(), tau0_host.stride_0(), w0_host.data());
783 Kokkos::deep_copy(A0_device, A0_host);
784 Kokkos::deep_copy(A0, A0_device);
785 Kokkos::deep_copy(tau0_device, tau0_host);
786 Kokkos::deep_copy(tau0, tau0_device);
788 Kokkos::parallel_for (systemName_,
789 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
790 KOKKOS_LAMBDA (
const size_t ic) {
791 auto w = Kokkos::subview(work, ic, Kokkos::ALL());
793 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
796 KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
797 1, A0.extent(0), A0.extent(1),
798 A0.data(), A0.stride_0(), A0.stride_1(),
799 tau0.data(), tau0.stride_0(),
800 b.data(), 1, b.stride_0(),
804 KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(
false,
807 A0.data(), A0.stride_0(), A0.stride_1(),
808 b.data(), b.stride_0());
811 for(ordinal_type i=0; i<n; ++i){
812 basisCoeffs(ic,elemDof(i)) = b(i);
818 Kokkos::parallel_for (systemName_,
819 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
820 KOKKOS_LAMBDA (
const size_t ic) {
822 auto A = Kokkos::subview(elemMat, ic, Kokkos::ALL(), Kokkos::ALL());
823 auto tau = Kokkos::subview(taul, ic, Kokkos::ALL());
824 auto w = Kokkos::subview(work, ic, Kokkos::ALL());
826 for(ordinal_type i=n; i<n+m; ++i)
827 for(ordinal_type j=0; j<n; ++j)
831 KokkosBatched::SerialQR_Internal::invoke(A.extent(0), A.extent(1),
832 A.data(), A.stride_0(), A.stride_1(), tau.data(), tau.stride_0(), w.data());
834 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
837 KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
838 1, A.extent(0), A.extent(1),
839 A.data(), A.stride_0(), A.stride_1(),
840 tau.data(), tau.stride_0(),
841 b.data(), 1, b.stride_0(),
845 KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(
false,
848 A.data(), A.stride_0(), A.stride_1(),
849 b.data(), b.stride_0());
852 for(ordinal_type i=0; i<n; ++i){
853 basisCoeffs(ic,elemDof(i)) = b(i);
863 template<
typename ViewType1,
typename ViewType2,
typename ViewType3,
typename ViewType4>
864 void solveHost(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 ,
865 ViewType3,
const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
866 using value_type =
typename ViewType2::value_type;
867 using device_type = DeviceType;
868 using host_exec_space = Kokkos::DefaultHostExecutionSpace;
869 using host_memory_space = Kokkos::HostSpace;
870 using host_device_type = Kokkos::Device<host_exec_space,host_memory_space>;
871 using vector_host_type = Kokkos::View<value_type*,host_device_type>;
872 using scratch_host_type = Kokkos::View<value_type*,host_exec_space::scratch_memory_space>;
873 using matrix_host_type = Kokkos::View<value_type**,Kokkos::LayoutLeft,host_device_type>;
874 using do_not_init_tag = Kokkos::ViewAllocateWithoutInitializing;
875 using host_team_policy_type = Kokkos::TeamPolicy<host_exec_space>;
876 using host_range_policy_type = Kokkos::RangePolicy<host_exec_space>;
882 const ordinal_type numCells = basisCoeffs.extent(0);
883 const ordinal_type numRows = m+n, numCols = n;
886 Teuchos::LAPACK<ordinal_type,value_type> lapack;
889 Kokkos::View<ordinal_type*,host_device_type> elemDof_host(do_not_init_tag(
"elemDof_host"), elemDof.extent(0));
891 auto elemDof_device = Kokkos::create_mirror_view(
typename device_type::memory_space(), elemDof_host);
892 Kokkos::deep_copy(elemDof_device, elemDof); Kokkos::fence();
893 Kokkos::deep_copy(elemDof_host, elemDof_device);
897 auto elemRhs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemRhs);
898 auto elemMat_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemMat);
901 auto basisCoeffs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), basisCoeffs);
903 if (matrixIndependentOfCell_) {
905 matrix_host_type A(do_not_init_tag(
"A"), numRows, numRows);
907 for (ordinal_type j=0;j<numRows;++j)
908 for (ordinal_type i=0;i<numRows;++i)
909 A(i, j) = elemMat_host(0, i, j);
911 for (ordinal_type j=0;j<numCols;++j)
912 for (ordinal_type i=numCols;i<numRows;++i)
916 ordinal_type lwork(-1);
918 ordinal_type info(0);
921 numRows, numRows, numCells,
922 nullptr, std::max(1,numRows),
923 nullptr, std::max(1,numRows),
929 matrix_host_type C(do_not_init_tag(
"C"), numRows, numCells);
931 host_range_policy_type policy(0, numCells);
934 (
"ProjectionTools::solveHost::matrixIndependentOfCell::pack",
935 policy, [=](
const ordinal_type & ic) {
936 for (ordinal_type i=0;i<numRows;++i)
937 C(i,ic) = elemRhs_host(ic, i);
942 vector_host_type work(do_not_init_tag(
"work"), lwork);
943 ordinal_type info(0);
945 numRows, numRows, numCells,
946 A.data(), std::max(1,numRows),
947 C.data(), std::max(1,numRows),
950 INTREPID2_TEST_FOR_EXCEPTION
951 (info != 0, std::runtime_error,
"GELS return non-zero info code");
955 (
"ProjectionTools::solveHost::matrixIndependentOfCell::unpack",
956 policy, [=](
const ordinal_type & ic) {
957 for (ordinal_type i=0;i<numCols;++i)
958 basisCoeffs_host(ic,elemDof_host(i)) = C(i,ic);
962 const ordinal_type level(0);
963 host_team_policy_type policy(numCells, 1, 1);
966 ordinal_type lwork(-1);
968 ordinal_type info(0);
972 nullptr, std::max(1,numRows),
973 nullptr, std::max(1,numRows),
979 const ordinal_type per_team_extent = numRows*numRows + numRows + lwork;
980 const ordinal_type per_team_scratch = scratch_host_type::shmem_size(per_team_extent);
981 policy.set_scratch_size(level, Kokkos::PerTeam(per_team_scratch));
985 (
"ProjectionTools::solveHost::matrixDependentOfCell",
986 policy, [=](
const typename host_team_policy_type::member_type& member) {
987 const ordinal_type ic = member.league_rank();
989 scratch_host_type scratch(member.team_scratch(level), per_team_extent);
990 value_type * sptr = scratch.data();
993 matrix_host_type A(sptr, numRows, numRows); sptr += A.span();
994 for (ordinal_type j=0;j<numRows;++j)
995 for (ordinal_type i=0;i<numRows;++i)
996 A(i, j) = elemMat_host(ic, i, j);
998 for (ordinal_type j=0;j<numCols;++j)
999 for (ordinal_type i=numCols;i<numRows;++i)
1002 vector_host_type c(sptr, numRows); sptr += c.span();
1003 for (ordinal_type i=0;i<numRows;++i)
1004 c(i) = elemRhs_host(ic, i);
1006 vector_host_type work(sptr, lwork); sptr += work.span();
1007 ordinal_type info(0);
1009 numRows, numRows, 1,
1010 A.data(), std::max(1,numRows),
1011 c.data(), std::max(1,numRows),
1014 INTREPID2_TEST_FOR_EXCEPTION
1015 (info != 0, std::runtime_error,
"GELS return non-zero info code");
1018 for (ordinal_type i=0;i<numCols;++i)
1019 basisCoeffs_host(ic,elemDof_host(i)) = c(i);
1022 Kokkos::deep_copy(basisCoeffs, basisCoeffs_host);
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_HEX_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_In_FEM class.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TRI_In_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_HEX_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_WEDGE_I1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_HEX_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
Header file for the Intrepid2::Experimental::ProjectionStruct.
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
An helper class to compute the evaluation points and weights needed for performing projections.
ordinal_type getNumTargetEvalPoints()
Returns number of points where to evaluate the target function.
void createL2ProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
Initialize the ProjectionStruct for L2 projections.