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  solveParallel(basisCoeffs, elemMat, elemRhs, tau,
652  w, elemDof, n, m);
653 #else
654  solveSerial(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 solveParallel(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 taul,
665  ViewType3 work,const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
666  using HostSpaceType = typename Kokkos::Impl::is_space<DeviceType>::host_mirror_space::execution_space;
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, HostSpaceType> 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, HostSpaceType> 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 solveSerial(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 ,
774  ViewType3, const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
775  using valueType = typename ViewType2::value_type;
776  using HostSpaceType = typename Kokkos::Impl::is_space<DeviceType>::host_mirror_space::execution_space;
777  Kokkos::View<valueType**,Kokkos::LayoutLeft,HostSpaceType>
778  serialElemMat("serialElemMat", n+m, n+m);
779  Teuchos::LAPACK<ordinal_type,valueType> lapack_;
780  ordinal_type numCells = basisCoeffs.extent(0);
781 
782  if(matrixIndependentOfCell_) {
783  ViewType2 elemRhsTrans("transRhs", elemRhs.extent(1), elemRhs.extent(0));
784  Kokkos::View<valueType**,Kokkos::LayoutLeft,HostSpaceType>
785  pivVec("pivVec", m+n + std::max(m+n, numCells), 1);
786 
787  Kokkos::View<valueType**,Kokkos::LayoutLeft,HostSpaceType> serialElemRhs("serialElemRhs", n+m, numCells);
788 
789  Kokkos::DynRankView<typename ViewType2::value_type, HostSpaceType> A_host("A0_host", elemMat.extent(1),elemMat.extent(2));
790  auto A_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), A_host);
791  Kokkos::deep_copy(A_device, Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL()));
792  Kokkos::deep_copy(A_host, A_device);
793 
794  auto b = Kokkos::create_mirror_view_and_copy(HostSpaceType(), elemRhs);
795 
796  auto serialBasisCoeffs = Kokkos::create_mirror_view_and_copy(
797  HostSpaceType(), basisCoeffs);
798 
799  for(ordinal_type i=0; i<m+n; ++i) {
800  for(ordinal_type ic=0; ic< numCells; ++ic)
801  serialElemRhs(i,ic) = b(ic,i);
802  for(ordinal_type j=0; j<n; ++j)
803  serialElemMat(j,i) = A_host(j,i);
804  }
805 
806  for(ordinal_type i=n; i<n+m; ++i)
807  for(ordinal_type j=0; j<n; ++j)
808  serialElemMat(i,j) = serialElemMat(j,i);
809 
810  ordinal_type info = 0;
811  lapack_.GELS('N', n+m, n+m, numCells,
812  serialElemMat.data(), serialElemMat.stride_1(),
813  serialElemRhs.data(), serialElemRhs.stride_1(),
814  pivVec.data(), pivVec.extent(0),
815  &info);
816 
817  for(ordinal_type i=0; i<n; ++i) {
818  for (ordinal_type ic = 0; ic < numCells; ic++)
819  serialBasisCoeffs(ic,elemDof(i)) = serialElemRhs(i,ic);
820  }
821  }
822  else {
823  Kokkos::View<valueType**,Kokkos::LayoutLeft,HostSpaceType> pivVec("pivVec", 2*(m+n), 1);
824  Kokkos::View<valueType**,Kokkos::LayoutLeft,HostSpaceType> serialElemRhs("serialElemRhs", n+m, 1 );
825  for (ordinal_type ic = 0; ic < numCells; ic++) {
826  auto A = Kokkos::create_mirror_view_and_copy(HostSpaceType(),
827  Kokkos::subview(elemMat, ic, Kokkos::ALL(), Kokkos::ALL()));
828  auto b = Kokkos::create_mirror_view_and_copy(HostSpaceType(),
829  Kokkos::subview(elemRhs, ic, Kokkos::ALL()));
830  auto basisCoeffs_ = Kokkos::subview(basisCoeffs, ic, Kokkos::ALL());
831  auto serialBasisCoeffs = Kokkos::create_mirror_view_and_copy(HostSpaceType(),
832  basisCoeffs_);
833 
834  Kokkos::deep_copy(serialElemMat,valueType(0)); //LAPACK might overwrite the matrix
835 
836  for(ordinal_type i=0; i<m+n; ++i) {
837  serialElemRhs(i,0) = b(i);
838  for(ordinal_type j=0; j<n; ++j)
839  serialElemMat(j,i) = A(j,i);
840  }
841 
842  for(ordinal_type i=n; i<n+m; ++i)
843  for(ordinal_type j=0; j<n; ++j)
844  serialElemMat(i,j) = serialElemMat(j,i);
845 
846  // Using GELS because the matrix can be close to singular.
847  ordinal_type info = 0;
848  lapack_.GELS('N', n+m, n+m, 1,
849  serialElemMat.data(), serialElemMat.stride_1(),
850  serialElemRhs.data(), serialElemRhs.stride_1(),
851  pivVec.data(), pivVec.extent(0),
852  &info);
853 
854  if (info) {
855  std::stringstream ss;
856  ss << ">>> ERROR (Intrepid::ProjectionTools::getBasisCoeffs): "
857  << "LAPACK return with error code: "
858  << info;
859  INTREPID2_TEST_FOR_EXCEPTION( true, std::runtime_error, ss.str().c_str() );
860  }
861 
862  for(ordinal_type i=0; i<n; ++i) {
863  serialBasisCoeffs(elemDof(i)) = serialElemRhs(i,0);
864  }
865  Kokkos::deep_copy(basisCoeffs_,serialBasisCoeffs);
866  }
867  }
868  }
869  };
870 
871 };
872 
873 } //Experimental
874 } //Intrepid2
875 
876 
877 // include templated function definitions
883 
884 #endif
885 
886 
887 
888 
889 
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.
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 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.
Class to solve a square system A x = b on each cell A is expected to be saddle a point (KKT) matrix o...
void solveSerial(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.