Compadre  1.5.9
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Compadre_LinearAlgebra_Declarations.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Compadre: COMpatible PArticle Discretization and REmap Toolkit
4 //
5 // Copyright 2018 NTESS and the Compadre contributors.
6 // SPDX-License-Identifier: BSD-2-Clause
7 // *****************************************************************************
8 // @HEADER
9 #ifndef _COMPADRE_LINEAR_ALGEBRA_DECLARATIONS_HPP_
10 #define _COMPADRE_LINEAR_ALGEBRA_DECLARATIONS_HPP_
11 
12 #include "Compadre_Config.h"
13 #include "Compadre_Typedefs.hpp"
15 
16 namespace Compadre {
17 
18 namespace GMLS_LinearAlgebra {
19 
20  /*! \brief Calculates two eigenvectors corresponding to two dominant eigenvalues
21  \param teamMember [in] - Kokkos::TeamPolicy member type (created by parallel_for)
22  \param V [out] - dimension * dimension Kokkos View
23  \param PtP [in] - dimension * dimension Kokkos View
24  \param dimensions [in] - dimension of PtP
25  */
26  KOKKOS_INLINE_FUNCTION
27  void largestTwoEigenvectorsThreeByThreeSymmetric(const member_type& teamMember, scratch_matrix_right_type V, scratch_matrix_right_type PtP, const int dimensions, pool_type& random_number_pool);
28 
29  /*! \brief Solves a batch of problems with QR+Pivoting
30 
31  ~ Note: Very strong assumption on B. ~
32 
33  A contains num_matrices * lda * nda data which is num_matrices different (lda x nda) matrices with valid entries of size (M x N), and
34  B contains num_matrices * ldb * ndb data which is num_matrices different (ldb x ndb) right hand sides.
35  B is assumed to have one of two forms:
36  - \f$M==N\f$, the valid entries are the first (N X NRHS)
37  - \f$M>N\f$, the valid entries are the first (NRHS)
38  i.e. for this case, B is intended to store non-zero entries from a diagonal matrix (as a vector). For the \f$k^{th}\f$ matrix, the (m,m) entry of a diagonal matrix would here be stored in the \f$m^{th}\f$ position.
39 
40 
41 
42  \param pm [in] - manager class for team and thread parallelism
43  \param A [in/out] - matrix A (in), meaningless workspace output (out)
44  \param lda [in] - row dimension of each matrix in A
45  \param nda [in] - columns dimension of each matrix in A
46  \param B [in/out] - right hand sides (in), solution (out)
47  \param ldb [in] - row dimension of each matrix in B
48  \param ndb [in] - column dimension of each matrix in B
49  \param M [in] - number of rows containing data (maximum rows over everything in batch) in A
50  \param N [in] - number of columns containing data in each matrix in A
51  \param NRHS [in] - number of columns containing data in each matrix in B
52  \param num_matrices [in] - number of problems
53  \param implicit_RHS [in] - determines whether RHS will be stored implicitly. If true, instead of RHS storing the full sqrt(W) explicitly, only the diagonal entries of sqrt(W) will be stored as a 1D array beginning at entry with matrix coordinate (0,0).
54  */
55  template <typename A_layout=layout_right, typename B_layout=layout_right, typename X_layout=layout_right>
56  void batchQRPivotingSolve(ParallelManager pm, double *A, int lda, int nda, double *B, int ldb, int ndb, int M, int N, int NRHS, const int num_matrices, const bool implicit_RHS = true);
57 
58 } // GMLS_LinearAlgebra
59 } // Compadre
60 
61 #endif
team_policy::member_type member_type
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.
void batchQRPivotingSolve(ParallelManager pm, double *A, int lda, int nda, double *B, int ldb, int ndb, int M, int N, int NRHS, const int num_matrices, const bool implicit_RHS)
Solves a batch of problems with QR+Pivoting.
Kokkos::View< double **, layout_right, Kokkos::MemoryTraits< Kokkos::Unmanaged > > scratch_matrix_right_type
Kokkos::Random_XorShift64_Pool pool_type