Intrepid2
Intrepid2_ProjectionTools.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_PROJECTIONTOOLS_HPP__
48#define __INTREPID2_PROJECTIONTOOLS_HPP__
49
50#include "Intrepid2_ConfigDefs.hpp"
51#include "Intrepid2_Types.hpp"
52#include "Intrepid2_Utils.hpp"
53
54#include "Shards_CellTopology.hpp"
55#include "Shards_BasicTopologies.hpp"
56
58
59#include "Intrepid2_Basis.hpp"
60
61// -- HGRAD family
65
68
69// -- HCURL family
72
76
77// -- HDIV family
80
84
85// -- Lower order family
88
91
95
99
100#include "Teuchos_LAPACK.hpp"
102
104
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"
110#endif
111
112namespace Intrepid2 {
113
114namespace Experimental {
115
116
117
181
182template<typename DeviceType>
184public:
185 using ExecSpaceType = typename DeviceType::execution_space;
186 using MemSpaceType = typename DeviceType::memory_space;
187 using EvalPointsType = typename ProjectionStruct<DeviceType, double>::EvalPointsType;
188
189
206 template<typename BasisType,
207 typename ortValueType, class ...ortProperties>
208 static void
209 getL2EvaluationPoints(typename BasisType::ScalarViewType evaluationPoints,
210 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
211 const BasisType* cellBasis,
213 const EvalPointsType evalPointType = EvalPointsType::TARGET
214 );
215
235 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
236 typename funValsValueType, class ...funValsProperties,
237 typename BasisType,
238 typename ortValueType, class ...ortProperties>
239 static void
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,
246
247
264 template<typename BasisType>
265 static void
266 getL2DGEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints,
267 const BasisType* cellBasis,
269 const EvalPointsType evalPointType = EvalPointsType::TARGET
270 );
271
295 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
296 typename funValsValueType, class ...funValsProperties,
297 typename BasisType,
298 typename ortValueType, class ...ortProperties>
299 static void
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,
305
328 template<typename basisViewType, typename targetViewType, typename BasisType>
329 static void
330 getL2DGBasisCoeffs(basisViewType basisCoeffs,
331 const targetViewType targetAtTargetEPoints,
332 const BasisType* cellBasis,
334
335
355 template<typename BasisType, typename OrientationViewType >
356 static void
357 getHGradEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints,
358 typename BasisType::ScalarViewType gradEvalPoints,
359 const OrientationViewType cellOrientations,
360 const BasisType* cellBasis,
362 const EvalPointsType evalPointType = EvalPointsType::TARGET
363 );
364
388 template<class BasisCoeffsViewType, class TargetValueViewType, class TargetGradViewType,
389 class BasisType, class OrientationViewType>
390 static void
391 getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs,
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,
399
400
420 template<typename BasisType,
421 typename ortValueType, class ...ortProperties>
422 static void
423 getHCurlEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints,
424 typename BasisType::ScalarViewType curlEvalPoints,
425 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
426 const BasisType* cellBasis,
428 const EvalPointsType evalPointType = EvalPointsType::TARGET
429 );
430
456 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
457 typename funValsValueType, class ...funValsProperties,
458 typename BasisType,
459 typename ortValueType, class ...ortProperties>
460 static void
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,
469
470
490 template<typename BasisType,
491 typename ortValueType, class ...ortProperties>
492 static void
493 getHDivEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints,
494 typename BasisType::ScalarViewType divEvalPoints,
495 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
496 const BasisType* cellBasis,
498 const EvalPointsType evalPointType = EvalPointsType::TARGET
499 );
500
524 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
525 typename funValsValueType, class ...funValsProperties,
526 typename BasisType,
527 typename ortValueType, class ...ortProperties>
528 static void
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,
537
554 template<typename BasisType,
555 typename ortValueType, class ...ortProperties>
556 static void
557 getHVolEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints,
558 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
559 const BasisType* cellBasis,
561 const EvalPointsType evalPointType = EvalPointsType::TARGET
562 );
563
582 template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
583 typename funValsValueType, class ...funValsProperties,
584 typename BasisType,
585 typename ortValueType, class ...ortProperties>
586 static void
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,
593
594
595
611 template<typename dstViewType,
612 typename dstBasisType,
613 typename srcViewType,
614 typename srcBasisType,
615 typename OrientationViewType>
616 static void
617 projectField(dstViewType dstCoeffs,
618 const dstBasisType* dstBasis,
619 const srcViewType srcCoeffs,
620 const srcBasisType* srcBasis,
621 const OrientationViewType cellOrientations){
622
623
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");
628
630 projStruct.createL2ProjectionStruct(dstBasis, srcBasis->getDegree());
631
632 ordinal_type numPoints = projStruct.getNumTargetEvalPoints();
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;
637
638 typename Kokkos::DynRankView<typename srcBasisType::PointValueType, DeviceType> evaluationPoints("evaluationPoints", numCells, numPoints, dim);
639 getL2EvaluationPoints(evaluationPoints,
640 cellOrientations,
641 dstBasis,
642 &projStruct);
643
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);
650 } else {
651 srcAtEvalPoints = outViewType("srcAtEvalPoints", numCells, numPoints);
652 refBasisAtEvalPoints = outViewType("refBasisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
653 basisAtEvalPoints = outViewType("basisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
654 }
655
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()));
658
659 // Modify basis values to account for orientations
661 refBasisAtEvalPoints,
662 cellOrientations,
663 srcBasis);
664
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);
671 }
672 }
673 });
674 ExecSpaceType().fence();
675
676 getL2BasisCoeffs(dstCoeffs,
677 srcAtEvalPoints,
678 evaluationPoints,
679 cellOrientations,
680 dstBasis,
681 &projStruct);
682 }
683
684
685
695 struct ElemSystem {
696
697
698 std::string systemName_;
699 bool matrixIndependentOfCell_;
700
707
708 ElemSystem (std::string systemName, bool matrixIndependentOfCell) :
709 systemName_(systemName), matrixIndependentOfCell_(matrixIndependentOfCell){};
710
711
712
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,
743 w, elemDof, n, m);
744#else
745 solveHost(basisCoeffs, elemMat, elemRhs, tau,
746 w, elemDof, n, m);
747#endif
748
749 }
750
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>;
758
759 ordinal_type numCells = basisCoeffs.extent(0);
760
761 if(matrixIndependentOfCell_) {
762 auto A0 = Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL());
763 auto tau0 = Kokkos::subview(taul, 0, Kokkos::ALL());
764
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);
769
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);
773
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()));
777
778 //computing QR of A0. QR is saved in A0 and tau0
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());
782
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);
787
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());
792
793 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
794
795 //b'*Q0 -> b
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(),
801 w.data());
802
803 // R0^{-1} b -> b
804 KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(false,
805 A0.extent(0),
806 1.0,
807 A0.data(), A0.stride_0(), A0.stride_1(),
808 b.data(), b.stride_0());
809
810 //scattering b into the basis coefficients
811 for(ordinal_type i=0; i<n; ++i){
812 basisCoeffs(ic,elemDof(i)) = b(i);
813 }
814 });
815
816 } else {
817
818 Kokkos::parallel_for (systemName_,
819 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
820 KOKKOS_LAMBDA (const size_t ic) {
821
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());
825
826 for(ordinal_type i=n; i<n+m; ++i)
827 for(ordinal_type j=0; j<n; ++j)
828 A(i,j) = A(j,i);
829
830 //computing QR of A. QR is saved in A and tau
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());
833
834 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
835
836 //b'*Q -> b
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(),
842 w.data());
843
844 // R^{-1} b -> b
845 KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(false,
846 A.extent(0),
847 1.0,
848 A.data(), A.stride_0(), A.stride_1(),
849 b.data(), b.stride_0());
850
851 //scattering b into the basis coefficients
852 for(ordinal_type i=0; i<n; ++i){
853 basisCoeffs(ic,elemDof(i)) = b(i);
854 }
855 });
856 }
857 }
858#endif
859
862
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>;
877
879 Kokkos::fence();
880
882 const ordinal_type numCells = basisCoeffs.extent(0);
883 const ordinal_type numRows = m+n, numCols = n;
884
886 Teuchos::LAPACK<ordinal_type,value_type> lapack;
887
889 Kokkos::View<ordinal_type*,host_device_type> elemDof_host(do_not_init_tag("elemDof_host"), elemDof.extent(0));
890 {
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);
894 }
895
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);
899
901 auto basisCoeffs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), basisCoeffs);
902
903 if (matrixIndependentOfCell_) {
905 matrix_host_type A(do_not_init_tag("A"), numRows, numRows);
906 {
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);
910
911 for (ordinal_type j=0;j<numCols;++j)
912 for (ordinal_type i=numCols;i<numRows;++i)
913 A(i, j) = A(j, i);
914 }
915
916 ordinal_type lwork(-1);
917 {
918 ordinal_type info(0);
919 value_type work[2];
920 lapack.GELS('N',
921 numRows, numRows, numCells,
922 nullptr, std::max(1,numRows),
923 nullptr, std::max(1,numRows),
924 &work[0], lwork,
925 &info);
926 lwork = work[0];
927 }
928
929 matrix_host_type C(do_not_init_tag("C"), numRows, numCells);
930
931 host_range_policy_type policy(0, numCells);
932 {
933 Kokkos::parallel_for
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);
938 });
939 }
940 {
942 vector_host_type work(do_not_init_tag("work"), lwork);
943 ordinal_type info(0);
944 lapack.GELS('N',
945 numRows, numRows, numCells,
946 A.data(), std::max(1,numRows),
947 C.data(), std::max(1,numRows),
948 work.data(), lwork,
949 &info);
950 INTREPID2_TEST_FOR_EXCEPTION
951 (info != 0, std::runtime_error, "GELS return non-zero info code");
952 }
953 {
954 Kokkos::parallel_for
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);
959 });
960 }
961 } else {
962 const ordinal_type level(0);
963 host_team_policy_type policy(numCells, 1, 1);
964
966 ordinal_type lwork(-1);
967 {
968 ordinal_type info(0);
969 value_type work[2];
970 lapack.GELS('N',
971 numRows, numRows, 1,
972 nullptr, std::max(1,numRows),
973 nullptr, std::max(1,numRows),
974 &work[0], lwork,
975 &info);
976 lwork = work[0];
977 }
978
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));
982
984 Kokkos::parallel_for
985 ("ProjectionTools::solveHost::matrixDependentOfCell",
986 policy, [=](const typename host_team_policy_type::member_type& member) {
987 const ordinal_type ic = member.league_rank();
988
989 scratch_host_type scratch(member.team_scratch(level), per_team_extent);
990 value_type * sptr = scratch.data();
991
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);
997
998 for (ordinal_type j=0;j<numCols;++j)
999 for (ordinal_type i=numCols;i<numRows;++i)
1000 A(i, j) = A(j, i);
1001
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);
1005
1006 vector_host_type work(sptr, lwork); sptr += work.span();
1007 ordinal_type info(0);
1008 lapack.GELS('N',
1009 numRows, numRows, 1,
1010 A.data(), std::max(1,numRows),
1011 c.data(), std::max(1,numRows),
1012 work.data(), lwork,
1013 &info);
1014 INTREPID2_TEST_FOR_EXCEPTION
1015 (info != 0, std::runtime_error, "GELS return non-zero info code");
1016
1018 for (ordinal_type i=0;i<numCols;++i)
1019 basisCoeffs_host(ic,elemDof_host(i)) = c(i);
1020 });
1021 }
1022 Kokkos::deep_copy(basisCoeffs, basisCoeffs_host);
1023 }
1024 };
1025
1026};
1027
1028} //Experimental
1029} //Intrepid2
1030
1031
1032// include templated function definitions
1038
1039#endif
1040
1041
1042
1043
1044
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::OrientationTools and Intrepid2::Impl::OrientationTools classes.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates,...
Header file for the Intrepid2::Experimental::ProjectionStruct.
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HCURL project...
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HDIV projecti...
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HGRAD project...
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HVOL projecti...
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for L2 projection...
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.
A class providing static members to perform projection-based interpolations:
static void getHVolBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, 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.
static void getL2BasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, 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 getL2DGEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for local L2 projection for broken HGRAD HCURL HDIV and HVOL spaces.
static void getHGradEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType gradEvalPoints, const OrientationViewType cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for HGrad projection.
static void getL2EvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for L2 projection.
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 getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs, const TargetValueViewType targetAtEvalPoints, const TargetGradViewType targetGradAtGradEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType gradEvalPoints, 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 getHCurlBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetCurlAtCurlEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType curlEvalPoints, 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 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 getHVolEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for HVol projection.
static void getHDivBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties... > basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties... > targetDivAtDivEvalPoints, const typename BasisType::ScalarViewType evaluationPoints, const typename BasisType::ScalarViewType divEvalPoints, 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 getHDivEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType divEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for HDiv projection.
static void getHCurlEvaluationPoints(typename BasisType::ScalarViewType evaluationPoints, typename BasisType::ScalarViewType curlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties... > cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct, const EvalPointsType evalPointType=EvalPointsType::TARGET)
Computes evaluation points for HCurl projection.
static void modifyBasisByOrientation(Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis)
Modify basis due to orientation.
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.
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...
ElemSystem(std::string systemName, bool matrixIndependentOfCell)
Functor constructor.