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 
57 #include "Intrepid2_PointTools.hpp"
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 
112 namespace Intrepid2 {
113 
114 namespace Experimental {
115 
116 
117 
182 template<typename DeviceType>
184 public:
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 
604  struct ElemSystem {
605 
606 
607  std::string systemName_;
608  bool matrixIndependentOfCell_;
609 
617  ElemSystem (std::string systemName, bool matrixIndependentOfCell) :
618  systemName_(systemName), matrixIndependentOfCell_(matrixIndependentOfCell){};
619 
620 
621 
647  template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
648  void solve(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 tau,
649  ViewType3 w,const ViewType4 elemDof, ordinal_type n, ordinal_type m=0) {
650 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
651  solveDevice(basisCoeffs, elemMat, elemRhs, tau,
652  w, elemDof, n, m);
653 #else
654  solveHost(basisCoeffs, elemMat, elemRhs, tau,
655  w, elemDof, n, m);
656 #endif
657 
658  }
659 
662 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
663  template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
664  void solveDevice(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 taul,
665  ViewType3 work,const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
666  using HostDeviceType = Kokkos::Device<Kokkos::DefaultHostExecutionSpace,Kokkos::HostSpace>;
667 
668  ordinal_type numCells = basisCoeffs.extent(0);
669 
670  if(matrixIndependentOfCell_) {
671  auto A0 = Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL());
672  auto tau0 = Kokkos::subview(taul, 0, Kokkos::ALL());
673 
674  Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> A0_host("A0_host", A0.extent(0),A0.extent(1));
675  auto A0_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), A0_host);
676  Kokkos::deep_copy(A0_device, A0);
677  Kokkos::deep_copy(A0_host, A0_device);
678 
679  for(ordinal_type i=n; i<n+m; ++i)
680  for(ordinal_type j=0; j<n; ++j)
681  A0_host(i,j) = A0_host(j,i);
682 
683  Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> tau0_host("A0_host", tau0.extent(0));
684  auto tau0_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), tau0_host);
685  auto w0_host = Kokkos::create_mirror_view(Kokkos::subview(work, 0, Kokkos::ALL()));
686 
687  //computing QR of A0. QR is saved in A0 and tau0
688  KokkosBatched::SerialQR_Internal::invoke(A0_host.extent(0), A0_host.extent(1),
689  A0_host.data(), A0_host.stride_0(), A0_host.stride_1(),
690  tau0_host.data(), tau0_host.stride_0(), w0_host.data());
691 
692  Kokkos::deep_copy(A0_device, A0_host);
693  Kokkos::deep_copy(A0, A0_device);
694  Kokkos::deep_copy(tau0_device, tau0_host);
695  Kokkos::deep_copy(tau0, tau0_device);
696 
697  Kokkos::parallel_for (systemName_,
698  Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
699  KOKKOS_LAMBDA (const size_t ic) {
700  auto w = Kokkos::subview(work, ic, Kokkos::ALL());
701 
702  auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
703 
704  //b'*Q0 -> b
705  KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
706  1, A0.extent(0), A0.extent(1),
707  A0.data(), A0.stride_0(), A0.stride_1(),
708  tau0.data(), tau0.stride_0(),
709  b.data(), 1, b.stride_0(),
710  w.data());
711 
712  // R0^{-1} b -> b
713  KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(false,
714  A0.extent(0),
715  1.0,
716  A0.data(), A0.stride_0(), A0.stride_1(),
717  b.data(), b.stride_0());
718 
719  //scattering b into the basis coefficients
720  for(ordinal_type i=0; i<n; ++i){
721  basisCoeffs(ic,elemDof(i)) = b(i);
722  }
723  });
724 
725  } else {
726 
727  Kokkos::parallel_for (systemName_,
728  Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
729  KOKKOS_LAMBDA (const size_t ic) {
730 
731  auto A = Kokkos::subview(elemMat, ic, Kokkos::ALL(), Kokkos::ALL());
732  auto tau = Kokkos::subview(taul, ic, Kokkos::ALL());
733  auto w = Kokkos::subview(work, ic, Kokkos::ALL());
734 
735  for(ordinal_type i=n; i<n+m; ++i)
736  for(ordinal_type j=0; j<n; ++j)
737  A(i,j) = A(j,i);
738 
739  //computing QR of A. QR is saved in A and tau
740  KokkosBatched::SerialQR_Internal::invoke(A.extent(0), A.extent(1),
741  A.data(), A.stride_0(), A.stride_1(), tau.data(), tau.stride_0(), w.data());
742 
743  auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
744 
745  //b'*Q -> b
746  KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
747  1, A.extent(0), A.extent(1),
748  A.data(), A.stride_0(), A.stride_1(),
749  tau.data(), tau.stride_0(),
750  b.data(), 1, b.stride_0(),
751  w.data());
752 
753  // R^{-1} b -> b
754  KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(false,
755  A.extent(0),
756  1.0,
757  A.data(), A.stride_0(), A.stride_1(),
758  b.data(), b.stride_0());
759 
760  //scattering b into the basis coefficients
761  for(ordinal_type i=0; i<n; ++i){
762  basisCoeffs(ic,elemDof(i)) = b(i);
763  }
764  });
765  }
766  }
767 #endif
768 
772  template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
773  void solveHost(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 ,
774  ViewType3, const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
775  using value_type = typename ViewType2::value_type;
776  using device_type = DeviceType;
777  using host_exec_space = Kokkos::DefaultHostExecutionSpace;
778  using host_memory_space = Kokkos::HostSpace;
779  using host_device_type = Kokkos::Device<host_exec_space,host_memory_space>;
780  using vector_host_type = Kokkos::View<value_type*,host_device_type>;
781  using scratch_host_type = Kokkos::View<value_type*,host_exec_space::scratch_memory_space>;
782  using matrix_host_type = Kokkos::View<value_type**,Kokkos::LayoutLeft,host_device_type>;
783  using do_not_init_tag = Kokkos::ViewAllocateWithoutInitializing;
784  using host_team_policy_type = Kokkos::TeamPolicy<host_exec_space>;
785  using host_range_policy_type = Kokkos::RangePolicy<host_exec_space>;
786 
788  Kokkos::fence();
789 
791  const ordinal_type numCells = basisCoeffs.extent(0);
792  const ordinal_type numRows = m+n, numCols = n;
793 
795  Teuchos::LAPACK<ordinal_type,value_type> lapack;
796 
798  Kokkos::View<ordinal_type*,host_device_type> elemDof_host(do_not_init_tag("elemDof_host"), elemDof.extent(0));
799  {
800  auto elemDof_device = Kokkos::create_mirror_view(typename device_type::memory_space(), elemDof_host);
801  Kokkos::deep_copy(elemDof_device, elemDof); Kokkos::fence();
802  Kokkos::deep_copy(elemDof_host, elemDof_device);
803  }
804 
806  auto elemRhs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemRhs);
807  auto elemMat_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemMat);
808 
810  auto basisCoeffs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), basisCoeffs);
811 
812  if (matrixIndependentOfCell_) {
814  matrix_host_type A(do_not_init_tag("A"), numRows, numRows);
815  {
816  for (ordinal_type j=0;j<numRows;++j)
817  for (ordinal_type i=0;i<numRows;++i)
818  A(i, j) = elemMat_host(0, i, j);
819 
820  for (ordinal_type j=0;j<numCols;++j)
821  for (ordinal_type i=numCols;i<numRows;++i)
822  A(i, j) = A(j, i);
823  }
824 
825  ordinal_type lwork(-1);
826  {
827  ordinal_type info(0);
828  value_type work[2];
829  lapack.GELS('N',
830  numRows, numRows, numCells,
831  nullptr, numRows,
832  nullptr, numRows,
833  &work[0], lwork,
834  &info);
835  lwork = work[0];
836  }
837 
838  matrix_host_type C(do_not_init_tag("C"), numRows, numCells);
839 
840  host_range_policy_type policy(0, numCells);
841  {
842  Kokkos::parallel_for
843  ("ProjectionTools::solveHost::matrixIndependentOfCell::pack",
844  policy, [=](const ordinal_type & ic) {
845  for (ordinal_type i=0;i<numRows;++i)
846  C(i,ic) = elemRhs_host(ic, i);
847  });
848  }
849  {
851  vector_host_type work(do_not_init_tag("work"), lwork);
852  ordinal_type info(0);
853  lapack.GELS('N',
854  numRows, numRows, numCells,
855  A.data(), numRows,
856  C.data(), numRows,
857  work.data(), lwork,
858  &info);
859  INTREPID2_TEST_FOR_EXCEPTION
860  (info != 0, std::runtime_error, "GELS return non-zero info code");
861  }
862  {
863  Kokkos::parallel_for
864  ("ProjectionTools::solveHost::matrixIndependentOfCell::unpack",
865  policy, [=](const ordinal_type & ic) {
866  for (ordinal_type i=0;i<numCols;++i)
867  basisCoeffs_host(ic,elemDof_host(i)) = C(i,ic);
868  });
869  }
870  } else {
871  const ordinal_type level(0);
872  host_team_policy_type policy(numCells, 1, 1);
873 
875  ordinal_type lwork(-1);
876  {
877  ordinal_type info(0);
878  value_type work[2];
879  lapack.GELS('N',
880  numRows, numRows, 1,
881  nullptr, numRows,
882  nullptr, numRows,
883  &work[0], lwork,
884  &info);
885  lwork = work[0];
886  }
887 
888  const ordinal_type per_team_extent = numRows*numRows + numRows + lwork;
889  const ordinal_type per_team_scratch = scratch_host_type::shmem_size(per_team_extent);
890  policy.set_scratch_size(level, Kokkos::PerTeam(per_team_scratch));
891 
893  Kokkos::parallel_for
894  ("ProjectionTools::solveHost::matrixDependentOfCell",
895  policy, [=](const typename host_team_policy_type::member_type& member) {
896  const ordinal_type ic = member.league_rank();
897 
898  scratch_host_type scratch(member.team_scratch(level), per_team_extent);
899  value_type * sptr = scratch.data();
900 
902  matrix_host_type A(sptr, numRows, numRows); sptr += A.span();
903  for (ordinal_type j=0;j<numRows;++j)
904  for (ordinal_type i=0;i<numRows;++i)
905  A(i, j) = elemMat_host(ic, i, j);
906 
907  for (ordinal_type j=0;j<numCols;++j)
908  for (ordinal_type i=numCols;i<numRows;++i)
909  A(i, j) = A(j, i);
910 
911  vector_host_type c(sptr, numRows); sptr += c.span();
912  for (ordinal_type i=0;i<numRows;++i)
913  c(i) = elemRhs_host(ic, i);
914 
915  vector_host_type work(sptr, lwork); sptr += work.span();
916  ordinal_type info(0);
917  lapack.GELS('N',
918  numRows, numRows, 1,
919  A.data(), numRows,
920  c.data(), numRows,
921  work.data(), lwork,
922  &info);
923  INTREPID2_TEST_FOR_EXCEPTION
924  (info != 0, std::runtime_error, "GELS return non-zero info code");
925 
927  for (ordinal_type i=0;i<numCols;++i)
928  basisCoeffs_host(ic,elemDof_host(i)) = c(i);
929  });
930  }
931  Kokkos::deep_copy(basisCoeffs, basisCoeffs_host);
932  }
933  };
934 
935 };
936 
937 } //Experimental
938 } //Intrepid2
939 
940 
941 // include templated function definitions
947 
948 #endif
949 
950 
951 
952 
953 
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.
Header file for the Intrepid2::Basis_HGRAD_LINE_Cn_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_HEX_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_WEDGE_I1_FEM class.
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 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.
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HCURL project...
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...
Header file for the Intrepid2::Basis_HDIV_HEX_In_FEM class.
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.
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for L2 projection...
Header function for Intrepid2::Util class and other utility functions.
Header file for the Intrepid2::Basis_HGRAD_TRI_Cn_FEM class.
Header file for the Intrepid2::Basis_HVOL_LINE_Cn_FEM class.
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...
Header file for the Intrepid2::Basis_HDIV_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_TET_In_FEM class.
A class providing static members to perform projection-based interpolations:
Header file for the Intrepid2::Basis_HGRAD_TET_Cn_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_QUAD_In_FEM class.
Header file for the Intrepid2::OrientationTools and Intrepid2::Impl::OrientationTools classes...
Header file for the Intrepid2::Basis_HCURL_QUAD_In_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_In_FEM class.
Class to solve a square system A x = b on each cell A is expected to be saddle a point (KKT) matrix o...
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HDIV projecti...
Contains definitions of custom data types in Intrepid2.
Header file for the Intrepid2::Basis_HGRAD_QUAD_Cn_FEM class.
Header file for the Intrepid2::Basis_HDIV_TET_In_FEM class.
Header file for the Intrepid2::Experimental::ProjectionStruct.
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HGRAD project...
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 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.
ElemSystem(std::string systemName, bool matrixIndependentOfCell)
Functor constructor.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_WEDGE_I1_FEM class.
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.
Header file for the Intrepid2::Basis_HCURL_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HCURL_HEX_I1_FEM class.
Header file for the Intrepid2::Experimental::ProjectionTools containing definitions for HVOL projecti...
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...
An helper class to compute the evaluation points and weights needed for performing projections...
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.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_FEM class.
Header file for the Intrepid2::Basis_HVOL_TRI_Cn_FEM class.
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 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.
Header file for the abstract base class Intrepid2::Basis.
Header file for the Intrepid2::Basis_HCURL_TRI_In_FEM class.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.
Header file for the Intrepid2::Basis_HGRAD_HEX_Cn_FEM class.