Compadre
1.2.0
|
Functions | |
void | batchQRFactorize (ParallelManager pm, double *P, int lda, int nda, double *RHS, int ldb, int ndb, int M, int N, int NRHS, const int num_matrices, const size_t max_neighbors=0, const int initial_index_of_batch=0, int *neighbor_list_sizes=NULL) |
Calls LAPACK or CUBLAS to solve a batch of QR problems. More... | |
void | batchSVDFactorize (ParallelManager pm, bool swap_layout_P, double *P, int lda, int nda, bool swap_layout_RHS, double *RHS, int ldb, int ndb, int M, int N, int NRHS, const int num_matrices, const size_t max_neighbors=0, const int initial_index_of_batch=0, int *neighbor_list_sizes=NULL) |
Calls LAPACK or CUBLAS to solve a batch of SVD problems. More... | |
void | batchLUFactorize (ParallelManager pm, double *P, int lda, int nda, double *RHS, int ldb, int ndb, int M, int N, int NRHS, const int num_matrices, const size_t max_neighbors=0, const int initial_index_of_batch=0, int *neighbor_list_sizes=NULL) |
Calls LAPACK or CUBLAS to solve a batch of LU problems. More... | |
KOKKOS_INLINE_FUNCTION void | createM (const member_type &teamMember, scratch_matrix_right_type M_data, scratch_matrix_right_type weighted_P, const int columns, const int rows) |
Creates a matrix M=A^T*A from a matrix A. More... | |
KOKKOS_INLINE_FUNCTION void | largestTwoEigenvectorsThreeByThreeSymmetric (const member_type &teamMember, scratch_matrix_right_type V, scratch_matrix_right_type PtP, const int dimensions, pool_type &random_number_pool) |
Calculates two eigenvectors corresponding to two dominant eigenvalues. More... | |
void Compadre::GMLS_LinearAlgebra::batchLUFactorize | ( | ParallelManager | pm, |
double * | P, | ||
int | lda, | ||
int | nda, | ||
double * | RHS, | ||
int | ldb, | ||
int | ndb, | ||
int | M, | ||
int | N, | ||
int | NRHS, | ||
const int | num_matrices, | ||
const size_t | max_neighbors = 0 , |
||
const int | initial_index_of_batch = 0 , |
||
int * | neighbor_list_sizes = NULL |
||
) |
Calls LAPACK or CUBLAS to solve a batch of LU problems.
P contains num_matrices * lda * ndb data which is num_matrices different matrices, and RHS contains num_matrices * ldn * ndb data which is num_matrices different matrix right hand sides.
pm | [in] - manager class for team and thread parallelism |
P | [in/out] - evaluation of sampling functional on polynomial basis (in), meaningless workspace output (out) |
lda | [in] - row dimension of each matrix in P |
nda | [in] - columns dimension of each matrix in P |
RHS | [in/out] - basis to invert P against (in), polynomial coefficients (out) |
ldb | [in] - row dimension of each matrix in RHS |
ndb | [in] - column dimension of each matrix in RHS |
M | [in] - number of rows containing data (maximum rows over everything in batch) |
N | [in] - number of columns containing data |
num_matrices | [in] - number of target (GMLS problems to be solved) |
max_neighbors | [in] - integer for maximum neighbor over all targets |
neighbor_list_sizes | [in] - pointer to all neighbor list sizes for each target |
Definition at line 581 of file Compadre_LinearAlgebra.cpp.
void Compadre::GMLS_LinearAlgebra::batchQRFactorize | ( | ParallelManager | pm, |
double * | P, | ||
int | lda, | ||
int | nda, | ||
double * | RHS, | ||
int | ldb, | ||
int | ndb, | ||
int | M, | ||
int | N, | ||
int | NRHS, | ||
const int | num_matrices, | ||
const size_t | max_neighbors = 0 , |
||
const int | initial_index_of_batch = 0 , |
||
int * | neighbor_list_sizes = NULL |
||
) |
Calls LAPACK or CUBLAS to solve a batch of QR problems.
P contains num_matrices * lda * ndb data which is num_matrices different matrices, and RHS contains num_matrices * ldn * ndb data which is num_matrices different matrix right hand sides.
pm | [in] - manager class for team and thread parallelism |
P | [in/out] - evaluation of sampling functional on polynomial basis (in), meaningless workspace output (out) |
lda | [in] - row dimension of each matrix in P |
nda | [in] - columns dimension of each matrix in P |
RHS | [in/out] - basis to invert P against (in), polynomial coefficients (out) |
ldb | [in] - row dimension of each matrix in RHS |
ndb | [in] - column dimension of each matrix in RHS |
M | [in] - number of rows containing data (maximum rows over everything in batch) |
N | [in] - number of columns containing data |
num_matrices | [in] - number of target (GMLS problems to be solved) |
max_neighbors | [in] - integer for maximum neighbor over all targets |
neighbor_list_sizes | [in] - pointer to all neighbor list sizes for each target |
Definition at line 14 of file Compadre_LinearAlgebra.cpp.
void Compadre::GMLS_LinearAlgebra::batchSVDFactorize | ( | ParallelManager | pm, |
bool | swap_layout_P, | ||
double * | P, | ||
int | lda, | ||
int | nda, | ||
bool | swap_layout_RHS, | ||
double * | RHS, | ||
int | ldb, | ||
int | ndb, | ||
int | M, | ||
int | N, | ||
int | NRHS, | ||
const int | num_matrices, | ||
const size_t | max_neighbors = 0 , |
||
const int | initial_index_of_batch = 0 , |
||
int * | neighbor_list_sizes = NULL |
||
) |
Calls LAPACK or CUBLAS to solve a batch of SVD problems.
P contains num_matrices * lda * ndb data which is num_matrices different matrices, and RHS contains num_matrices * ldn * ndb data which is num_matrices different matrix right hand sides.
pm | [in] - manager class for team and thread parallelism |
swap_layout_P | [in] - boolean to check if layout of P will need to be swapped or not |
P | [in/out] - evaluation of sampling functional on polynomial basis (in), meaningless workspace output (out) |
lda | [in] - row dimension of each matrix in P |
nda | [in] - columns dimension of each matrix in P |
swap_layout_RHS | [in] - boolean to check if layout of RHS will need to be swapped or not |
RHS | [in/out] - basis to invert P against (in), polynomial coefficients (out) |
ldb | [in] - row dimension of each matrix in RHS |
ndb | [in] - column dimension of each matrix in RHS |
M | [in] - number of rows containing data (maximum rows over everything in batch) |
N | [in] - number of columns containing data |
num_matrices | [in] - number of target (GMLS problems to be solved) |
max_neighbors | [in] - integer for maximum neighbor over all targets |
neighbor_list_sizes | [in] - pointer to all neighbor list sizes for each target |
Definition at line 164 of file Compadre_LinearAlgebra.cpp.
KOKKOS_INLINE_FUNCTION void Compadre::GMLS_LinearAlgebra::createM | ( | const member_type & | teamMember, |
scratch_matrix_right_type | M_data, | ||
scratch_matrix_right_type | weighted_P, | ||
const int | columns, | ||
const int | rows | ||
) |
Creates a matrix M=A^T*A from a matrix A.
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
M_data | [out] - result of weighted_P^T * weighted_P |
weighted_P | [in] - matrix to be multiplied by its transpose |
columns | [in] - number of columns in weighted_P |
rows | [in] - number of rows in weighted_P |
Definition at line 10 of file Compadre_LinearAlgebra_Definitions.hpp.
KOKKOS_INLINE_FUNCTION void Compadre::GMLS_LinearAlgebra::largestTwoEigenvectorsThreeByThreeSymmetric | ( | const member_type & | teamMember, |
scratch_matrix_right_type | V, | ||
scratch_matrix_right_type | PtP, | ||
const int | dimensions, | ||
pool_type & | random_number_pool | ||
) |
Calculates two eigenvectors corresponding to two dominant eigenvalues.
teamMember | [in] - Kokkos::TeamPolicy member type (created by parallel_for) |
V | [out] - dimension * dimension Kokkos View |
PtP | [in] - dimension * dimension Kokkos View |
dimensions | [in] - dimension of PtP |
Definition at line 55 of file Compadre_LinearAlgebra_Definitions.hpp.