Intrepid2
Public Member Functions | Public Attributes | List of all members
Intrepid2::ProjectionTools< DeviceType >::ElemSystem Struct Reference

Class to solve a square system A x = b on each cell A is expected to be saddle a point (KKT) matrix of the form [C B; B^T 0], where C has size nxn and B nxm, with n>0, m>=0. B^T is copied from B, so one does not have to define the B^T portion of A. b will contain the solution x. The first n-entries of x are copied into the provided basis coefficients using the provided indexing. The system is solved either with a QR factorization implemented in KokkosKernels or with Lapack GELS function. More...

#include <Intrepid2_ProjectionTools.hpp>

Public Member Functions

 ElemSystem (std::string systemName, bool matrixIndependentOfCell)
 Functor constructor. More...
 
template<typename ViewType1 , typename ViewType2 , typename ViewType3 , typename ViewType4 >
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 Lapack GELS depending on whether Kokkos Kernel is enabled. More...
 
template<typename ViewType1 , typename ViewType2 , typename ViewType3 , typename ViewType4 >
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. More...
 

Public Attributes

std::string systemName_
 
bool matrixIndependentOfCell_
 

Detailed Description

template<typename DeviceType>
struct Intrepid2::ProjectionTools< DeviceType >::ElemSystem

Class to solve a square system A x = b on each cell A is expected to be saddle a point (KKT) matrix of the form [C B; B^T 0], where C has size nxn and B nxm, with n>0, m>=0. B^T is copied from B, so one does not have to define the B^T portion of A. b will contain the solution x. The first n-entries of x are copied into the provided basis coefficients using the provided indexing. The system is solved either with a QR factorization implemented in KokkosKernels or with Lapack GELS function.

Definition at line 472 of file Intrepid2_ProjectionTools.hpp.

Constructor & Destructor Documentation

template<typename DeviceType >
Intrepid2::ProjectionTools< DeviceType >::ElemSystem::ElemSystem ( std::string  systemName,
bool  matrixIndependentOfCell 
)
inline

Functor constructor.

Parameters
systemName[in] - string containing the name of the system (passed to parallel for)
matrixIndependentOfCell[in] - bool: whether the local cell matrix of the system changes from cell to cell if true, the matrix factorization is preformed only on the first cell and reused on other cells.

Definition at line 485 of file Intrepid2_ProjectionTools.hpp.

Member Function Documentation

template<typename DeviceType >
template<typename ViewType1 , typename ViewType2 , typename ViewType3 , typename ViewType4 >
void Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve ( ViewType1  basisCoeffs,
ViewType2  elemMat,
ViewType2  elemRhs,
ViewType2  tau,
ViewType3  w,
const ViewType4  elemDof,
ordinal_type  n,
ordinal_type  m = 0 
)
inline

Solve the system and returns the basis coefficients solve the system either using Kokkos Kernel QR or Lapack GELS depending on whether Kokkos Kernel is enabled.

C - num. cells
P - num. evaluation points
Parameters
basisCoeffs[out] - rank-2 view (C,F) containing the basis coefficients
elemMat[in/out] - rank-3 view (C,P,P) containing the element matrix of size numCells x (n+m)x(n+m) on each cell it will be overwritten.
elemRhs[in/out] - rank-2 view (C,P) containing the element rhs on each cell of size numCells x (n+m) it will contain the solution of the system on output
tau[out] - rank-2 view (C,P) used to store the QR factorization size: numCells x (n+m)
w[out] - rank-2 view (C,P) used has a workspace Layout Right, size: numCells x (n+m)
elemDof[in] - rank-1 view having dimension n, containing the basis numbering
n[in] - ordinal_type, basis cardinality
m[in] - ordinal_type, dimension of the constraint of the KKT system

Definition at line 516 of file Intrepid2_ProjectionTools.hpp.

References Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solveHost().

Referenced by Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHGradBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs(), and Intrepid2::ProjectionTools< DeviceType >::getL2DGBasisCoeffs().

template<typename DeviceType >
template<typename ViewType1 , typename ViewType2 , typename ViewType3 , typename ViewType4 >
void Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solveHost ( ViewType1  basisCoeffs,
ViewType2  elemMat,
ViewType2  elemRhs,
ViewType2  ,
ViewType3  ,
const ViewType4  elemDof,
ordinal_type  n,
ordinal_type  m 
)
inline

Parallel implementation of solve, using Kokkos Kernels QR factoriation.

Serial implementation of solve, using Lapack GELS function

make sure all on-going kernels are done

const values

capture without this pointer

stride view copy

mirror view to host

this in-out variable

invert the first matrix and apply for all

workspace query

GELS does scaling and separating QR and apply Q is not stable

workspace query

solve for all

either A comes from host or device it is not column major; it needs repack

scatter back to system

Definition at line 641 of file Intrepid2_ProjectionTools.hpp.

Referenced by Intrepid2::ProjectionTools< DeviceType >::ElemSystem::solve().


The documentation for this struct was generated from the following file: