14 #ifndef __INTREPID2_PROJECTIONTOOLS_HPP__
15 #define __INTREPID2_PROJECTIONTOOLS_HPP__
17 #include "Intrepid2_ConfigDefs.hpp"
21 #include "Shards_CellTopology.hpp"
22 #include "Shards_BasicTopologies.hpp"
45 #include "Teuchos_LAPACK.hpp"
50 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
51 #include "KokkosBatched_QR_Serial_Internal.hpp"
52 #include "KokkosBatched_ApplyQ_Serial_Internal.hpp"
53 #include "KokkosBatched_Trsv_Serial_Internal.hpp"
54 #include "KokkosBatched_Util.hpp"
111 template<
typename DeviceType>
114 using ExecSpaceType =
typename DeviceType::execution_space;
115 using MemSpaceType =
typename DeviceType::memory_space;
116 using EvalPointsType =
typename ProjectionStruct<DeviceType, double>::EvalPointsType;
138 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
139 typename funValsValueType,
class ...funValsProperties,
141 typename ortValueType,
class ...ortProperties>
143 getL2BasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
144 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
145 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
146 const BasisType* cellBasis,
173 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
174 typename funValsValueType,
class ...funValsProperties,
176 typename ortValueType,
class ...ortProperties>
178 getL2DGBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
179 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
180 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
181 const BasisType* cellBasis,
207 template<
typename basisViewType,
typename targetViewType,
typename BasisType>
210 const targetViewType targetAtTargetEPoints,
211 const BasisType* cellBasis,
235 template<
class BasisCoeffsViewType,
class TargetValueViewType,
class TargetGradViewType,
236 class BasisType,
class OrientationViewType>
239 const TargetValueViewType targetAtEvalPoints,
240 const TargetGradViewType targetGradAtGradEvalPoints,
241 const OrientationViewType cellOrientations,
242 const BasisType* cellBasis,
269 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
270 typename funValsValueType,
class ...funValsProperties,
272 typename ortValueType,
class ...ortProperties>
274 getHCurlBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
275 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
276 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetCurlAtCurlEvalPoints,
277 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
278 const BasisType* cellBasis,
302 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
303 typename funValsValueType,
class ...funValsProperties,
305 typename ortValueType,
class ...ortProperties>
307 getHDivBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
308 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
309 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetDivAtDivEvalPoints,
310 const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
311 const BasisType* cellBasis,
332 template<
typename basisCoeffsValueType,
class ...basisCoeffsProperties,
333 typename funValsValueType,
class ...funValsProperties,
335 typename ortValueType,
class ...ortProperties>
337 getHVolBasisCoeffs(Kokkos::DynRankView<basisCoeffsValueType,basisCoeffsProperties...> basisCoeffs,
338 const Kokkos::DynRankView<funValsValueType,funValsProperties...> targetAtEvalPoints,
339 [[maybe_unused]]
const Kokkos::DynRankView<ortValueType, ortProperties...> cellOrientations,
340 const BasisType* cellBasis,
360 template<
typename dstViewType,
361 typename dstBasisType,
362 typename srcViewType,
363 typename srcBasisType,
364 typename OrientationViewType>
367 const dstBasisType* dstBasis,
368 const srcViewType srcCoeffs,
369 const srcBasisType* srcBasis,
370 const OrientationViewType cellOrientations){
373 INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getFunctionSpace() != srcBasis->getFunctionSpace(), std::runtime_error,
374 "The source and destination bases are not compatible. They need to belong to the same function space");
375 INTREPID2_TEST_FOR_EXCEPTION(dstBasis->getBaseCellTopology().getKey() != srcBasis->getBaseCellTopology().getKey(), std::runtime_error,
376 "The source and destination bases are not compatible. They do not have the same basic cell topology");
382 ordinal_type numCells = cellOrientations.extent(0);
383 ordinal_type dim = srcBasis->getBaseCellTopology().getDimension();
384 ordinal_type srcBasisCardinality = srcBasis->getCardinality();
385 ordinal_type fieldDimension = (srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HCURL || srcBasis->getFunctionSpace() == Intrepid2::FUNCTION_SPACE_HDIV) ? dim : 1;
387 auto evaluationPoints = projStruct.getAllEvalPoints();
388 ordinal_type numPoints = evaluationPoints.extent(0);
390 using outViewType = Kokkos::DynRankView<typename srcBasisType::OutputValueType, DeviceType>;
391 outViewType srcAtEvalPoints, refBasisAtEvalPoints, basisAtEvalPoints;
392 if(fieldDimension == dim) {
393 srcAtEvalPoints = outViewType(
"srcAtEvalPoints", numCells, numPoints, dim);
394 refBasisAtEvalPoints = outViewType(
"refBasisAtEvalPoints", srcBasisCardinality, numPoints, dim);
395 basisAtEvalPoints = outViewType(
"basisAtEvalPoints", numCells, srcBasisCardinality, numPoints, dim);
397 srcAtEvalPoints = outViewType(
"srcAtEvalPoints", numCells, numPoints);
398 refBasisAtEvalPoints = outViewType(
"refBasisAtEvalPoints", srcBasisCardinality, numPoints);
399 basisAtEvalPoints = outViewType(
"basisAtEvalPoints", numCells, srcBasisCardinality, numPoints);
402 srcBasis->getValues(refBasisAtEvalPoints,evaluationPoints);
406 refBasisAtEvalPoints,
410 Kokkos::parallel_for(Kokkos::RangePolicy<typename DeviceType::execution_space>(0,numCells),
411 KOKKOS_LAMBDA (
const int &ic) {
412 for(
int j=0; j<numPoints; ++j) {
413 for(
int k=0; k<srcBasisCardinality; ++k) {
414 for(
int d=0; d<fieldDimension; ++d)
415 srcAtEvalPoints.access(ic,j,d) += srcCoeffs(ic,k)*basisAtEvalPoints.access(ic,k,j,d);
419 ExecSpaceType().fence();
442 std::string systemName_;
443 bool matrixIndependentOfCell_;
452 ElemSystem (std::string systemName,
bool matrixIndependentOfCell) :
453 systemName_(systemName), matrixIndependentOfCell_(matrixIndependentOfCell){};
482 template<
typename ViewType1,
typename ViewType2,
typename ViewType3,
typename ViewType4>
483 void solve(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 tau,
484 ViewType3 w,
const ViewType4 elemDof, ordinal_type n, ordinal_type m=0) {
485 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
486 solveDevice(basisCoeffs, elemMat, elemRhs, tau,
489 solveHost(basisCoeffs, elemMat, elemRhs, tau,
497 #ifdef HAVE_INTREPID2_KOKKOSKERNELS
498 template<
typename ViewType1,
typename ViewType2,
typename ViewType3,
typename ViewType4>
499 void solveDevice(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 taul,
500 ViewType3 work,
const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
501 using HostDeviceType = Kokkos::Device<Kokkos::DefaultHostExecutionSpace,Kokkos::HostSpace>;
503 ordinal_type numCells = basisCoeffs.extent(0);
505 if(matrixIndependentOfCell_) {
506 auto A0 = Kokkos::subview(elemMat, 0, Kokkos::ALL(), Kokkos::ALL());
507 auto tau0 = Kokkos::subview(taul, 0, Kokkos::ALL());
509 Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> A0_host(
"A0_host", A0.extent(0),A0.extent(1));
510 auto A0_device = Kokkos::create_mirror_view(
typename DeviceType::memory_space(), A0_host);
511 Kokkos::deep_copy(A0_device, A0);
512 Kokkos::deep_copy(A0_host, A0_device);
514 for(ordinal_type i=n; i<n+m; ++i)
515 for(ordinal_type j=0; j<n; ++j)
516 A0_host(i,j) = A0_host(j,i);
518 Kokkos::DynRankView<typename ViewType2::value_type, HostDeviceType> tau0_host(
"A0_host", tau0.extent(0));
519 auto tau0_device = Kokkos::create_mirror_view(
typename DeviceType::memory_space(), tau0_host);
520 auto w0_host = Kokkos::create_mirror_view(Kokkos::subview(work, 0, Kokkos::ALL()));
523 KokkosBatched::SerialQR_Internal::invoke(A0_host.extent(0), A0_host.extent(1),
524 A0_host.data(), A0_host.stride_0(), A0_host.stride_1(),
525 tau0_host.data(), tau0_host.stride_0(), w0_host.data());
527 Kokkos::deep_copy(A0_device, A0_host);
528 Kokkos::deep_copy(A0, A0_device);
529 Kokkos::deep_copy(tau0_device, tau0_host);
530 Kokkos::deep_copy(tau0, tau0_device);
532 Kokkos::parallel_for (systemName_,
533 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
534 KOKKOS_LAMBDA (
const size_t ic) {
535 auto w = Kokkos::subview(work, ic, Kokkos::ALL());
537 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
540 KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
541 1, A0.extent(0), A0.extent(1),
542 A0.data(), A0.stride_0(), A0.stride_1(),
543 tau0.data(), tau0.stride_0(),
544 b.data(), 1, b.stride_0(),
548 KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(
false,
551 A0.data(), A0.stride_0(), A0.stride_1(),
552 b.data(), b.stride_0());
555 for(ordinal_type i=0; i<n; ++i){
556 basisCoeffs(ic,elemDof(i)) = b(i);
562 Kokkos::parallel_for (systemName_,
563 Kokkos::RangePolicy<ExecSpaceType, int> (0, numCells),
564 KOKKOS_LAMBDA (
const size_t ic) {
566 auto A = Kokkos::subview(elemMat, ic, Kokkos::ALL(), Kokkos::ALL());
567 auto tau = Kokkos::subview(taul, ic, Kokkos::ALL());
568 auto w = Kokkos::subview(work, ic, Kokkos::ALL());
570 for(ordinal_type i=n; i<n+m; ++i)
571 for(ordinal_type j=0; j<n; ++j)
575 KokkosBatched::SerialQR_Internal::invoke(A.extent(0), A.extent(1),
576 A.data(), A.stride_0(), A.stride_1(), tau.data(), tau.stride_0(), w.data());
578 auto b = Kokkos::subview(elemRhs, ic, Kokkos::ALL());
581 KokkosBatched::SerialApplyQ_RightForwardInternal::invoke(
582 1, A.extent(0), A.extent(1),
583 A.data(), A.stride_0(), A.stride_1(),
584 tau.data(), tau.stride_0(),
585 b.data(), 1, b.stride_0(),
589 KokkosBatched::SerialTrsvInternalUpper<KokkosBatched::Algo::Trsv::Unblocked>::invoke(
false,
592 A.data(), A.stride_0(), A.stride_1(),
593 b.data(), b.stride_0());
596 for(ordinal_type i=0; i<n; ++i){
597 basisCoeffs(ic,elemDof(i)) = b(i);
607 template<
typename ViewType1,
typename ViewType2,
typename ViewType3,
typename ViewType4>
608 void solveHost(ViewType1 basisCoeffs, ViewType2 elemMat, ViewType2 elemRhs, ViewType2 ,
609 ViewType3,
const ViewType4 elemDof, ordinal_type n, ordinal_type m) {
610 using value_type =
typename ViewType2::value_type;
611 using device_type = DeviceType;
612 using host_exec_space = Kokkos::DefaultHostExecutionSpace;
613 using host_memory_space = Kokkos::HostSpace;
614 using host_device_type = Kokkos::Device<host_exec_space,host_memory_space>;
615 using vector_host_type = Kokkos::View<value_type*,host_device_type>;
616 using scratch_host_type = Kokkos::View<value_type*,host_exec_space::scratch_memory_space>;
617 using matrix_host_type = Kokkos::View<value_type**,Kokkos::LayoutLeft,host_device_type>;
618 using do_not_init_tag = Kokkos::ViewAllocateWithoutInitializing;
619 using host_team_policy_type = Kokkos::TeamPolicy<host_exec_space>;
620 using host_range_policy_type = Kokkos::RangePolicy<host_exec_space>;
626 const ordinal_type numCells = basisCoeffs.extent(0);
627 const ordinal_type numRows = m+n, numCols = n;
630 Teuchos::LAPACK<ordinal_type,value_type> lapack;
633 Kokkos::View<ordinal_type*,host_device_type> elemDof_host(do_not_init_tag(
"elemDof_host"), elemDof.extent(0));
635 auto elemDof_device = Kokkos::create_mirror_view(
typename device_type::memory_space(), elemDof_host);
636 Kokkos::deep_copy(elemDof_device, elemDof); Kokkos::fence();
637 Kokkos::deep_copy(elemDof_host, elemDof_device);
641 auto elemRhs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemRhs);
642 auto elemMat_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), elemMat);
645 auto basisCoeffs_host = Kokkos::create_mirror_view_and_copy(host_memory_space(), basisCoeffs);
647 if (matrixIndependentOfCell_) {
649 matrix_host_type A(do_not_init_tag(
"A"), numRows, numRows);
651 for (ordinal_type j=0;j<numRows;++j)
652 for (ordinal_type i=0;i<numRows;++i)
653 A(i, j) = elemMat_host(0, i, j);
655 for (ordinal_type j=0;j<numCols;++j)
656 for (ordinal_type i=numCols;i<numRows;++i)
660 ordinal_type lwork(-1);
662 ordinal_type info(0);
665 numRows, numRows, numCells,
666 nullptr, std::max(1,numRows),
667 nullptr, std::max(1,numRows),
673 matrix_host_type C(do_not_init_tag(
"C"), numRows, numCells);
675 host_range_policy_type policy(0, numCells);
678 (
"ProjectionTools::solveHost::matrixIndependentOfCell::pack",
679 policy, [=](
const ordinal_type & ic) {
680 for (ordinal_type i=0;i<numRows;++i)
681 C(i,ic) = elemRhs_host(ic, i);
686 vector_host_type work(do_not_init_tag(
"work"), lwork);
687 ordinal_type info(0);
689 numRows, numRows, numCells,
690 A.data(), std::max(1,numRows),
691 C.data(), std::max(1,numRows),
694 INTREPID2_TEST_FOR_EXCEPTION
695 (info != 0, std::runtime_error,
"GELS return non-zero info code");
699 (
"ProjectionTools::solveHost::matrixIndependentOfCell::unpack",
700 policy, [=](
const ordinal_type & ic) {
701 for (ordinal_type i=0;i<numCols;++i)
702 basisCoeffs_host(ic,elemDof_host(i)) = C(i,ic);
706 const ordinal_type level(0);
707 host_team_policy_type policy(numCells, 1, 1);
710 ordinal_type lwork(-1);
712 ordinal_type info(0);
716 nullptr, std::max(1,numRows),
717 nullptr, std::max(1,numRows),
723 const ordinal_type per_team_extent = numRows*numRows + numRows + lwork;
724 const ordinal_type per_team_scratch = scratch_host_type::shmem_size(per_team_extent);
725 policy.set_scratch_size(level, Kokkos::PerTeam(per_team_scratch));
729 (
"ProjectionTools::solveHost::matrixDependentOfCell",
730 policy, [=](
const typename host_team_policy_type::member_type& member) {
731 const ordinal_type ic = member.league_rank();
733 scratch_host_type scratch(member.team_scratch(level), per_team_extent);
734 value_type * sptr = scratch.data();
737 matrix_host_type A(sptr, numRows, numRows); sptr += A.span();
738 for (ordinal_type j=0;j<numRows;++j)
739 for (ordinal_type i=0;i<numRows;++i)
740 A(i, j) = elemMat_host(ic, i, j);
742 for (ordinal_type j=0;j<numCols;++j)
743 for (ordinal_type i=numCols;i<numRows;++i)
746 vector_host_type c(sptr, numRows); sptr += c.span();
747 for (ordinal_type i=0;i<numRows;++i)
748 c(i) = elemRhs_host(ic, i);
750 vector_host_type work(sptr, lwork); sptr += work.span();
751 ordinal_type info(0);
754 A.data(), std::max(1,numRows),
755 c.data(), std::max(1,numRows),
758 INTREPID2_TEST_FOR_EXCEPTION
759 (info != 0, std::runtime_error,
"GELS return non-zero info code");
762 for (ordinal_type i=0;i<numCols;++i)
763 basisCoeffs_host(ic,elemDof_host(i)) = c(i);
766 Kokkos::deep_copy(basisCoeffs, basisCoeffs_host);
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 function for Intrepid2::Util class and other utility functions.
Header file for the Intrepid2::Basis_HDIV_QUAD_I1_FEM class.
Header file for the Intrepid2::Basis_HDIV_TRI_I1_FEM class.
Contains definitions of custom data types in Intrepid2.
Header file for the Intrepid2::ProjectionStruct.
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.
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::Basis_HCURL_TRI_I1_FEM class.
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...
Header file for the abstract base class Intrepid2::Basis.