Intrepid2
Intrepid2_ProjectionTools.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
14 #ifndef __INTREPID2_PROJECTIONTOOLS_HPP__
15 #define __INTREPID2_PROJECTIONTOOLS_HPP__
16 
17 #include "Intrepid2_ConfigDefs.hpp"
18 #include "Intrepid2_Types.hpp"
19 #include "Intrepid2_Utils.hpp"
20 
21 #include "Shards_CellTopology.hpp"
22 #include "Shards_BasicTopologies.hpp"
23 
24 #include "Intrepid2_PointTools.hpp"
25 
26 #include "Intrepid2_Basis.hpp"
27 
29 
30 // -- Lower order family
33 
36 
40 
44 
45 #include "Teuchos_LAPACK.hpp"
47 
49 
50 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
51 #include "KokkosBatched_QR_Serial_Internal.hpp"
52 #include "KokkosBatched_ApplyQ_Serial_Internal.hpp"
53 #if KOKKOS_VERSION >= 40599
54 #include "KokkosBatched_Trsv_Decl.hpp"
55 #else
56 #include "KokkosBatched_Trsv_Serial_Internal.hpp"
57 #endif
58 #include "KokkosBatched_Util.hpp"
59 #endif
60 
61 namespace Intrepid2 {
62 
115 template<typename DeviceType>
117 public:
118  using ExecSpaceType = typename DeviceType::execution_space;
119  using MemSpaceType = typename DeviceType::memory_space;
120  using EvalPointsType = typename ProjectionStruct<DeviceType, double>::EvalPointsType;
121 
122 
142  template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
143  typename funValsValueType, class ...funValsProperties,
144  typename BasisType,
145  typename ortValueType, class ...ortProperties>
146  static void
147  getL2BasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
148  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
149  const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
150  const BasisType* cellBasis,
152 
153 
177  template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
178  typename funValsValueType, class ...funValsProperties,
179  typename BasisType,
180  typename ortValueType, class ...ortProperties>
181  static void
182  getL2DGBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
183  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
184  const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
185  const BasisType* cellBasis,
187 
188 
211  template<typename basisViewType, typename targetViewType, typename BasisType>
212  static void
213  getL2DGBasisCoeffs(basisViewType basisCoeffs,
214  const targetViewType targetAtTargetEPoints,
215  const BasisType* cellBasis,
217 
239  template<class BasisCoeffsViewType, class TargetValueViewType, class TargetGradViewType,
240  class BasisType, class OrientationViewType>
241  static void
242  getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs,
243  const TargetValueViewType targetAtEvalPoints,
244  const TargetGradViewType targetGradAtGradEvalPoints,
245  const OrientationViewType cellOrientations,
246  const BasisType* cellBasis,
248 
249 
273  template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
274  typename funValsValueType, class ...funValsProperties,
275  typename BasisType,
276  typename ortValueType, class ...ortProperties>
277  static void
278  getHCurlBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
279  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
280  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetCurlAtCurlEvalPoints,
281  const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
282  const BasisType* cellBasis,
284 
306  template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
307  typename funValsValueType, class ...funValsProperties,
308  typename BasisType,
309  typename ortValueType, class ...ortProperties>
310  static void
311  getHDivBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
312  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
313  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetDivAtDivEvalPoints,
314  const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
315  const BasisType* cellBasis,
317 
318 
336  template<typename basisCoeffsValueType, class ...basisCoeffsProperties,
337  typename funValsValueType, class ...funValsProperties,
338  typename BasisType,
339  typename ortValueType, class ...ortProperties>
340  static void
341  getHVolBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
342  const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
343  [[maybe_unused]] const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
344  const BasisType* cellBasis,
346 
347 
348 
364  template<typename dstViewType,
365  typename dstBasisType,
366  typename srcViewType,
367  typename srcBasisType,
368  typename OrientationViewType>
369  static void
370  projectField(dstViewType dstCoeffs,
371  const dstBasisType* dstBasis,
372  const srcViewType srcCoeffs,
373  const srcBasisType* srcBasis,
374  const OrientationViewType cellOrientations){
375 
376 
377  INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getFunctionSpace() != srcBasis->getFunctionSpace(), std::runtime_error,
378  "The source and destination bases are not compatible. They need to belong to the same function space");
379  INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getBaseCellTopology().getKey() != srcBasis->getBaseCellTopology().getKey(), std::runtime_error,
380  "The source and destination bases are not compatible. They do not have the same basic cell topology");
381 
383  projStruct.createL2ProjectionStruct(dstBasis, srcBasis->getDegree());
384 
385 
386  ordinal_type numCells = cellOrientations.extent(0);
387  ordinal_type dim = srcBasis->getBaseCellTopology().getDimension();
388  ordinal_type srcBasisCardinality = srcBasis->getCardinality();
389  ordinal_type fieldDimension = (srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HCURL || srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HDIV) ? dim : 1;
390 
391  auto evaluationPoints = projStruct.getAllEvalPoints();
392  ordinal_type numPoints = evaluationPoints.extent(0);
393 
394  using outViewType = Kokkos::DynRankView<typename srcBasisType::OutputValueType, DeviceType>;
395  outViewType srcAtEvalPoints, refBasisAtEvalPoints, basisAtEvalPoints;
396  if(fieldDimension == dim) {
397  srcAtEvalPoints = outViewType("srcAtEvalPoints", numCells, numPoints, dim);
398  refBasisAtEvalPoints = outViewType("refBasisAtEvalPoints", srcBasisCardinality, numPoints, dim);
399  basisAtEvalPoints = outViewType("basisAtEvalPoints", numCells, srcBasisCardinality, numPoints, dim);
400  } else {
401  srcAtEvalPoints = outViewType("srcAtEvalPoints", numCells, numPoints);
402  refBasisAtEvalPoints = outViewType("refBasisAtEvalPoints", srcBasisCardinality, numPoints);
403  basisAtEvalPoints = outViewType("basisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
404  }
405 
406  srcBasis->getValues(refBasisAtEvalPoints,evaluationPoints);
407 
408  // Modify basis values to account for orientations
410  refBasisAtEvalPoints,
411  cellOrientations,
412  srcBasis);
413 
414  Kokkos::parallel_for(Kokkos::RangePolicy<typename DeviceType::execution_space>(0,numCells),
415  KOKKOS_LAMBDA (const int &ic) {
416  for(int j=0; j<numPoints; ++j) {
417  for(int k=0; k<srcBasisCardinality; ++k) {
418  for(int d=0; d<fieldDimension; ++d)
419  srcAtEvalPoints.access(ic,j,d) += srcCoeffs(ic,k)*basisAtEvalPoints.access(ic,k,j,d);
420  }
421  }
422  });
423  ExecSpaceType().fence();
424 
425  getL2BasisCoeffs(dstCoeffs,
426  srcAtEvalPoints,
427  cellOrientations,
428  dstBasis,
429  &projStruct);
430  }
431 
432 
433 
443  struct ElemSystem {
444 
445 
446  std::string systemName_;
447  bool matrixIndependentOfCell_;
448 
456  ElemSystem (std::string systemName, bool matrixIndependentOfCell) :
457  systemName_(systemName), matrixIndependentOfCell_(matrixIndependentOfCell){};
458 
459 
460 
486  template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
487  void solve(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 tau,
488  ViewType3 w,const ViewType4 elemDof, ordinal_type n, ordinal_type m=0) {
489 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
490  solveDevice(basisCoeffs, elemMat, elemRhs, tau,
491  w, elemDof, n, m);
492 #else
493  solveHost(basisCoeffs, elemMat, elemRhs, tau,
494  w, elemDof, n, m);
495 #endif
496 
497  }
498 
501 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
502  template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
503  void solveDevice(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 taul,
504  ViewType3 work,const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
505  using HostDeviceType = Kokkos::Device<Kokkos::DefaultHostExecutionSpace,Kokkos::HostSpace>;
506 
507  ordinal_type numCells = basisCoeffs.extent(0);
508 
509  if(matrixIndependentOfCell_) {
510  auto A0 = Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL());
511  auto tau0 = Kokkos::subview(taul, 0, Kokkos::ALL());
512 
513  Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> A0_host("A0_host", A0.extent(0),A0.extent(1));
514  auto A0_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), A0_host);
515  Kokkos::deep_copy(A0_device, A0);
516  Kokkos::deep_copy(A0_host, A0_device);
517 
518  for(ordinal_type i=n; i<n+m; ++i)
519  for(ordinal_type j=0; j<n; ++j)
520  A0_host(i,j) = A0_host(j,i);
521 
522  Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> tau0_host("A0_host", tau0.extent(0));
523  auto tau0_device = Kokkos::create_mirror_view(typename DeviceType::memory_space(), tau0_host);
524  auto w0_host = Kokkos::create_mirror_view(Kokkos::subview(work, 0, Kokkos::ALL()));
525 
526  //computing QR of A0. QR is saved in A0 and tau0
527  KokkosBatched::SerialQR_Internal::invoke(A0_host.extent(0), A0_host.extent(1),
528  A0_host.data(), A0_host.stride_0(), A0_host.stride_1(),
529  tau0_host.data(), tau0_host.stride_0(), w0_host.data());
530 
531  Kokkos::deep_copy(A0_device, A0_host);
532  Kokkos::deep_copy(A0, A0_device);
533  Kokkos::deep_copy(tau0_device, tau0_host);
534  Kokkos::deep_copy(tau0, tau0_device);
535 
536  Kokkos::parallel_for (systemName_,
537  Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
538  KOKKOS_LAMBDA (const size_t ic) {
539  auto w = Kokkos::subview(work, ic, Kokkos::ALL());
540 
541  auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
542 
543  //b'*Q0 -> b
544  KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
545  1, A0.extent(0), A0.extent(1),
546  A0.data(), A0.stride_0(), A0.stride_1(),
547  tau0.data(), tau0.stride_0(),
548  b.data(), 1, b.stride_0(),
549  w.data());
550 
551  // R0^{-1} b -> b
552 #if KOKKOS_VERSION >= 40599
553  KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(1.0, A0, b);
554 #else
555  KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(false,
556  A0.extent(0),
557  1.0,
558  A0.data(), A0.stride_0(), A0.stride_1(),
559  b.data(), b.stride_0());
560 #endif
561 
562  //scattering b into the basis coefficients
563  for(ordinal_type i=0; i<n; ++i){
564  basisCoeffs(ic,elemDof(i)) = b(i);
565  }
566  });
567 
568  } else {
569 
570  Kokkos::parallel_for (systemName_,
571  Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
572  KOKKOS_LAMBDA (const size_t ic) {
573 
574  auto A = Kokkos::subview(elemMat, ic, Kokkos::ALL(), Kokkos::ALL());
575  auto tau = Kokkos::subview(taul, ic, Kokkos::ALL());
576  auto w = Kokkos::subview(work, ic, Kokkos::ALL());
577 
578  for(ordinal_type i=n; i<n+m; ++i)
579  for(ordinal_type j=0; j<n; ++j)
580  A(i,j) = A(j,i);
581 
582  //computing QR of A. QR is saved in A and tau
583  KokkosBatched::SerialQR_Internal::invoke(A.extent(0), A.extent(1),
584  A.data(), A.stride_0(), A.stride_1(), tau.data(), tau.stride_0(), w.data());
585 
586  auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
587 
588  //b'*Q -> b
589  KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
590  1, A.extent(0), A.extent(1),
591  A.data(), A.stride_0(), A.stride_1(),
592  tau.data(), tau.stride_0(),
593  b.data(), 1, b.stride_0(),
594  w.data());
595 
596  // R^{-1} b -> b
597 #if KOKKOS_VERSION >= 40599
598  KokkosBatched::SerialTrsv<KokkosBatched::Uplo::Upper, KokkosBatched::Trans::NoTranspose, KokkosBatched::Diag::NonUnit, KokkosBatched::Algo::Trsv::Unblocked>::invoke(1.0, A, b);
599 #else
600  KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(false,
601  A.extent(0),
602  1.0,
603  A.data(), A.stride_0(), A.stride_1(),
604  b.data(), b.stride_0());
605 #endif
606 
607  //scattering b into the basis coefficients
608  for(ordinal_type i=0; i<n; ++i){
609  basisCoeffs(ic,elemDof(i)) = b(i);
610  }
611  });
612  }
613  }
614 #endif
615 
619  template<typename ViewType1, typename ViewType2, typename ViewType3, typename ViewType4>
620  void solveHost(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 ,
621  ViewType3, const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
622  using value_type = typename ViewType2::value_type;
623  using device_type = DeviceType;
624  using host_exec_space = Kokkos::DefaultHostExecutionSpace;
625  using host_memory_space = Kokkos::HostSpace;
626  using host_device_type = Kokkos::Device<host_exec_space,host_memory_space>;
627  using vector_host_type = Kokkos::View<value_type*,host_device_type>;
628  using scratch_host_type = Kokkos::View<value_type*,host_exec_space::scratch_memory_space>;
629  using matrix_host_type = Kokkos::View<value_type**,Kokkos::LayoutLeft,host_device_type>;
630  using do_not_init_tag = Kokkos::ViewAllocateWithoutInitializing;
631  using host_team_policy_type = Kokkos::TeamPolicy<host_exec_space>;
632  using host_range_policy_type = Kokkos::RangePolicy<host_exec_space>;
633 
635  Kokkos::fence();
636 
638  const ordinal_type numCells = basisCoeffs.extent(0);
639  const ordinal_type numRows = m+n, numCols = n;
640 
642  Teuchos::LAPACK<ordinal_type,value_type> lapack;
643 
645  Kokkos::View<ordinal_type*,host_device_type> elemDof_host(do_not_init_tag("elemDof_host"), elemDof.extent(0));
646  {
647  auto elemDof_device = Kokkos::create_mirror_view(typename device_type::memory_space(), elemDof_host);
648  Kokkos::deep_copy(elemDof_device, elemDof); Kokkos::fence();
649  Kokkos::deep_copy(elemDof_host, elemDof_device);
650  }
651 
653  auto elemRhs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemRhs);
654  auto elemMat_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemMat);
655 
657  auto basisCoeffs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), basisCoeffs);
658 
659  if (matrixIndependentOfCell_) {
661  matrix_host_type A(do_not_init_tag("A"), numRows, numRows);
662  {
663  for (ordinal_type j=0;j<numRows;++j)
664  for (ordinal_type i=0;i<numRows;++i)
665  A(i, j) = elemMat_host(0, i, j);
666 
667  for (ordinal_type j=0;j<numCols;++j)
668  for (ordinal_type i=numCols;i<numRows;++i)
669  A(i, j) = A(j, i);
670  }
671 
672  ordinal_type lwork(-1);
673  {
674  ordinal_type info(0);
675  value_type work[2];
676  lapack.GELS('N',
677  numRows, numRows, numCells,
678  nullptr, std::max(1,numRows),
679  nullptr, std::max(1,numRows),
680  &work[0], lwork,
681  &info);
682  lwork = work[0];
683  }
684 
685  matrix_host_type C(do_not_init_tag("C"), numRows, numCells);
686 
687  host_range_policy_type policy(0, numCells);
688  {
689  Kokkos::parallel_for
690  ("ProjectionTools::solveHost::matrixIndependentOfCell::pack",
691  policy, [=](const ordinal_type & ic) {
692  for (ordinal_type i=0;i<numRows;++i)
693  C(i,ic) = elemRhs_host(ic, i);
694  });
695  }
696  {
698  vector_host_type work(do_not_init_tag("work"), lwork);
699  ordinal_type info(0);
700  lapack.GELS('N',
701  numRows, numRows, numCells,
702  A.data(), std::max(1,numRows),
703  C.data(), std::max(1,numRows),
704  work.data(), lwork,
705  &info);
706  INTREPID2_TEST_FOR_EXCEPTION
707  (info != 0, std::runtime_error, "GELS return non-zero info code");
708  }
709  {
710  Kokkos::parallel_for
711  ("ProjectionTools::solveHost::matrixIndependentOfCell::unpack",
712  policy, [=](const ordinal_type & ic) {
713  for (ordinal_type i=0;i<numCols;++i)
714  basisCoeffs_host(ic,elemDof_host(i)) = C(i,ic);
715  });
716  }
717  } else {
718  const ordinal_type level(0);
719  host_team_policy_type policy(numCells, 1, 1);
720 
722  ordinal_type lwork(-1);
723  {
724  ordinal_type info(0);
725  value_type work[2];
726  lapack.GELS('N',
727  numRows, numRows, 1,
728  nullptr, std::max(1,numRows),
729  nullptr, std::max(1,numRows),
730  &work[0], lwork,
731  &info);
732  lwork = work[0];
733  }
734 
735  const ordinal_type per_team_extent = numRows*numRows + numRows + lwork;
736  const ordinal_type per_team_scratch = scratch_host_type::shmem_size(per_team_extent);
737  policy.set_scratch_size(level, Kokkos::PerTeam(per_team_scratch));
738 
740  Kokkos::parallel_for
741  ("ProjectionTools::solveHost::matrixDependentOfCell",
742  policy, [=](const typename host_team_policy_type::member_type& member) {
743  const ordinal_type ic = member.league_rank();
744 
745  scratch_host_type scratch(member.team_scratch(level), per_team_extent);
746  value_type * sptr = scratch.data();
747 
749  matrix_host_type A(sptr, numRows, numRows); sptr += A.span();
750  for (ordinal_type j=0;j<numRows;++j)
751  for (ordinal_type i=0;i<numRows;++i)
752  A(i, j) = elemMat_host(ic, i, j);
753 
754  for (ordinal_type j=0;j<numCols;++j)
755  for (ordinal_type i=numCols;i<numRows;++i)
756  A(i, j) = A(j, i);
757 
758  vector_host_type c(sptr, numRows); sptr += c.span();
759  for (ordinal_type i=0;i<numRows;++i)
760  c(i) = elemRhs_host(ic, i);
761 
762  vector_host_type work(sptr, lwork); sptr += work.span();
763  ordinal_type info(0);
764  lapack.GELS('N',
765  numRows, numRows, 1,
766  A.data(), std::max(1,numRows),
767  c.data(), std::max(1,numRows),
768  work.data(), lwork,
769  &info);
770  INTREPID2_TEST_FOR_EXCEPTION
771  (info != 0, std::runtime_error, "GELS return non-zero info code");
772 
774  for (ordinal_type i=0;i<numCols;++i)
775  basisCoeffs_host(ic,elemDof_host(i)) = c(i);
776  });
777  }
778  Kokkos::deep_copy(basisCoeffs, basisCoeffs_host);
779  }
780  };
781 
782 };
783 
784 } //Intrepid2
785 
786 
787 // include templated function definitions
793 
794 #endif
795 
796 
797 
798 
799 
static void getHCurlBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetCurlAtCurlEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HCurl projection of the target function.
static void getHGradBasisCoeffs(BasisCoeffsViewType basisCoeffs, const TargetValueViewType targetAtEvalPoints, const TargetGradViewType targetGradAtGradEvalPoints, const OrientationViewType cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HGrad projection of the target function.
Header file for the Intrepid2::Basis_HDIV_TET_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_HEX_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_WEDGE_I1_FEM class.
Header file for the Intrepid2::ProjectionTools containing definitions for HCURL projections.
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...
Header file for the Intrepid2::ProjectionTools containing definitions for L2 projections.
static void getHDivBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetDivAtDivEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the HDiv projection of the target function.
Header function for Intrepid2::Util class and other utility functions.
Header file for the Intrepid2::Basis_HDIV_QUAD_I1_FEM class.
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.
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...
A class providing static members to perform projection-based interpolations:
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
static void modifyBasisByOrientation(Kokkos::DynRankView< outputValueType, outputProperties...> output, const Kokkos::DynRankView< inputValueType, inputProperties...> input, const OrientationViewType orts, const BasisType *basis, const bool transpose=false)
Modify basis due to orientation.
Header file for the Intrepid2::OrientationTools and Intrepid2::Impl::OrientationTools classes...
Header file for the Intrepid2::ProjectionTools containing definitions for HDIV projections.
Contains definitions of custom data types in Intrepid2.
ElemSystem(std::string systemName, bool matrixIndependentOfCell)
Functor constructor.
Header file for the Intrepid2::ProjectionStruct.
Header file for the Intrepid2::ProjectionTools containing definitions for HGRAD projections.
void createL2ProjectionStruct(const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
Initialize the ProjectionStruct for L2 projections.
Header file for the Intrepid2::Basis_HCURL_QUAD_I1_FEM class.
static void getHVolBasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, [[maybe_unused]] 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.
Header file for the Intrepid2::Basis_HCURL_WEDGE_I1_FEM class.
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::ProjectionTools containing definitions for HVOL projections.
static void getL2BasisCoeffs(Kokkos::DynRankView< basisCoeffsValueType, basisCoeffsProperties...> basisCoeffs, const Kokkos::DynRankView< funValsValueType, funValsProperties...> targetAtEvalPoints, const Kokkos::DynRankView< ortValueType, ortProperties...> cellOrientations, const BasisType *cellBasis, ProjectionStruct< DeviceType, typename BasisType::scalarType > *projStruct)
Computes the basis coefficients of the L2 projection of the target function.
Header file for the Intrepid2::Basis_HCURL_TRI_I1_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...
An helper class to compute the evaluation points and weights needed for performing projections...
Stateless class that acts as a factory for a family of nodal bases (hypercube topologies only at this...
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...
Header file for the abstract base class Intrepid2::Basis.
Header file for Intrepid2::PointTools class to provide utilities for barycentric coordinates, equispaced lattices, and warp-blend point distrubtions.