ChaseMpiProperties¶
-
template<class
T
>
classchase::mpi
::
ChaseMpiProperties
¶ A class to setup MPI properties and matrix distribution scheme for the implementation of ChASE on distributed-memory systems.
The ChaseMpiProperties class creates a 2D grid of MPI nodes with fixed width and height. It defines also the scheme to distribute an Hermitian matrix
A
of sideN
over this 2D grid of MPI nodes. Currently, two data distribution fashion are supported in ChASE.Block distribution scheme in which one submatrix of
A
is assigned to one single MPI node;Block-Cyclic distribution scheme, which distributes a series of submatrices of
A
to the MPI nodes in a round-robin manner so that each MPI rank gets several non-adjacent blocks. More details about this distribution can be found on Netlib.- Template Parameters
T
: the scalar type used for the application. ChASE is templated for real and complex scalars with both Single Precision and Double Precision, thusT
can be one offloat
,double
,std::complex<float>
andstd::complex<double>
.
Public Functions
-
ChaseMpiProperties
(std::size_t N, std::size_t mb, std::size_t nb, std::size_t nev, std::size_t nex, int row_dim, int col_dim, char *grid_major, int irsrc, int icsrc, MPI_Comm comm)¶ A constructor of the class ChaseMpiProperties which distributes matrix
A
in Block-Cyclic Distribution.It constructs a 2D grid of MPI ranks within the MPI communicator
comm_
.The dimension of this 2D grid is determined by explicit values
row_dim
andcol_dim
, thusrow_dim * col_dim
.The major to index the coordinations of the MPI ranks is determined by the value of
grid_major
: ifgrid_major="C"
, they are indexed bycolumn major
, ifgrid_major="R"
, they are indexed byRow major
.It distributes the Hermitian matrix
A
in a Block-Dsitribution scheme.
This constructor requires the explicit values for the initalization of the size
N
of the matrix A, the number of sought after extremal eigenvaluesnev
, and the number of extra eigenvaluenex
which defines, together withnev
, the search space, and the working MPI communicatorcomm_
for the construction of eigenproblem. For the distribution in block-cyclic scheme, it requires also the dimensions of 2D gridrow_dim
andcol_dim
, the row and column blocking factorsmb
andnb
, the major of grid indexinggrid_major
, and the process row/column over which the first row/column of the global matrixA
is distributed.The parameters related to Block-Cyclic Distribution conforms equally with the Array Descriptor of in-core Dense Matrices in ScaLAPACK.
All the private members are either initialized directly by these parameters, or setup within the construction of this constructor.
- Parameters
N
: Size of the square matrix defining the eigenproblem.mb
: Row blocking factor for Block-Cyclic Distribution.nb
: Column blocking factor for Block-Cyclic Distribution`.nev
: Number of desired extremal eigenvalues.nex
: Number of eigenvalues augmenting the search space. Usually a relatively small fraction ofnev
.row_dim
: number of row in the 2D grid of MPI nodes to be constructed.col_dim
: number of column in the 2D grid of MPI nodes to be constructed.grid_major
: the type of numbering of the MPI ranks within the constructed 2D grid.Row major numbering if
grid_major=="C"
.Column major numbering if
grid_major=="R"
.
irsrc
: Process row over which the first row of the global matrixA
is distributed.icsrc
: Process column over which the first column of the global matrixA
is distributed.comm
: the working MPI communicator for ChASE.
-
ChaseMpiProperties
(std::size_t N, std::size_t nev, std::size_t nex, std::size_t m, std::size_t n, int npr, int npc, char *grid_major, MPI_Comm comm)¶ A constructor of the class ChaseMpiProperties which distributes matrix
A
inBlock Distribution
.It constructs a 2D grid of MPI ranks within the MPI communicator
comm_
.The dimensions of this 2D grid is determined by the input arguments
npr
andnpc
. The 2D grid isnpr x npc
It distributes the Hermitian matrix
A
in a Block-Dsitribution scheme.
This constructor requires the explicit values for the initalization of the size
N
of the matrix A, the number of sought after extremal eigenvaluesnev
, and the number of extra eigenvaluenex
which defines, together withnev
, the search space, the dimension of local matrixm
andn
, the 2D MPI gridnpr
andnpc
, and the working MPI communicatorcomm_
.All the private members are either initialized directly by these parameters, or setup within the construction of this constructor.
- Parameters
N
: Size of the square matrix defining the eigenproblem.nev
: Number of desired extremal eigenvalues.nex
: Number of eigenvalues augmenting the search space. Usually a relatively small fraction ofnev
.m
: row number of local matrix on each MPI rankn
: column number of local matrix on each MPI ranknpr
: row number of 2D MPI gridnpc
: column number of 2D MPI gridgrid_major
: the type of numbering of the MPI ranks within the constructed 2D grid.Row major numbering if
grid_major=="C"
.Column major numbering if
grid_major=="R"
.
comm
: the working MPI communicator for ChASE.
-
ChaseMpiProperties
(std::size_t N, std::size_t nev, std::size_t nex, MPI_Comm comm)¶ A constructor of the class ChaseMpiProperties which distributes matrix
A
inBlock Distribution
.It constructs a 2D grid of MPI ranks within the MPI communicator
comm_
.The dimensions of this 2D grid is determined by
MPI_Cart_create
using all the available MPI nodes withincomm_
.It distributes the Hermitian matrix
A
in a Block-Dsitribution scheme.
This constructor requires the explicit values for the initalization of the size
N
of the matrix A, the number of sought after extremal eigenvaluesnev
, and the number of extra eigenvaluenex
which defines, together withnev
, the search space, and the working MPI communicatorcomm_
.All the private members are either initialized directly by these parameters, or setup within the construction of this constructor.
- Parameters
N
: Size of the square matrix defining the eigenproblem.nev
: Number of desired extremal eigenvalues.nex
: Number of eigenvalues augmenting the search space. Usually a relatively small fraction ofnev
.comm
: the working MPI communicator for ChASE.
-
std::size_t
get_N
()¶ Returns the rank of matrix
A
which is distributed within 2D MPI grid.- Return
N_
: the rank of matrixA
.
-
std::size_t
get_ldh
()¶ Return leading dimension of local part of Symmetric/Hermtian matrix on each MPI proc.
- Return
ldh_
, a private member of this class.
-
std::size_t
get_n
()¶ Returns column number of the local matrix on each MPI node.
- Return
n_
: the column number of the local matrix on each MPI node.
-
std::size_t
get_m
()¶ Returns row number of the local matrix.
- Return
m_
: Row number of the local matrix on each MPI node.
-
std::size_t
get_max_block
()¶ Returns the maximum column number of rectangular matrix
V
.- Return
max_block_
: Maximum column number of matrixV
.
-
std::size_t
GetNev
()¶ Returns number of desired extremal eigenpairs, which was set by users.
- Return
nev_
: Number of desired extremal eigenpairs.
-
std::size_t
GetNex
()¶ Returns the Increment of the search subspace so that its total size, which was set by users within chase::ChaseConfig.
- Return
nex_
: Increment of the search subspace.
-
std::size_t
get_nb
()¶ Returns the column blocking factor, this function is useful for ChASE with Block-Cyclic Distribution.
- Return
nb_
: Column blocking factor.
-
std::size_t
get_mb
()¶ Returns the row blocking factor, this function is useful for ChASE with Block-Cyclic Distribution.
- Return
mb_
: Row blocking factor.
-
MPI_Comm
get_row_comm
()¶ - Return
row_comm_
: the row communicator within 2D MPI grid.
-
MPI_Comm
get_col_comm
()¶ - Return
col_comm_
: the column communicator within 2D MPI grid.
-
int *
get_dims
()¶ Returns the dimension of cartesian communicator grid.
- Return
dims_
: an array of size 2 which store the dimensions of the constructed 2D MPI grid.
-
std::size_t *
get_off
()¶ This member function is mostly used for ChASE with Block-Distribution.
Returns offset of row and column of the local matrix on each MPI node regarding the index of row and column of global matrix
A
.For example, for a global matrix
A
,A_sub
is a submatrix ofA
of sizem * n
, which is extracted out ofA
starting from the row indexingi
and column indexingj
. This member function returnsi
andj
as an array of size 2 for the submatrix on each MPI rank.- Return
off_
: an array of size 2 which stores the offset of row and column of the local matrix on each MPI node regarding the index of row and column of global matrixA
.
-
void
get_off
(std::size_t *xoff, std::size_t *yoff, std::size_t *xlen, std::size_t *ylen)¶ This member function is mostly used for ChASE with Block-Distribution.
Return the offset and length along the two dimensions of the local matrix on each MPI node regarding the global index of matrix
A
.For example, for a global matrix
A
,A_sub
is a submatrix ofA
of sizem * n
, which is extracted out ofA
starting from the row indexingi
and column indexingj
. This member function outputsxoff=i
,yoff=j
,xlen=m
, andylen=n
for the submatrix on each MPI rank.- Parameters
[in/out]
:xoff
-> the offset of row of the local matrix on each MPI node regarding the row index of matrixA
.[in/out]
:yoff
-> the offset of column of the local matrix on each MPI node regarding the column index of matrixA
.[in/out]
:xlen
-> the number of row of the local matrix on each MPI node.[in/out]
:ylen
-> the number of column of the local matrix on each MPI node.
-
std::size_t
get_mblocks
()¶ This member function is mostly used for ChASE with Block-Cyclic Distribution. The number of submatrix on each MPI node is
mblocks_ * nblocks_
.- Return
mblocks_
: the number of submatrices in the local matrix on each MPI node along the row direction.
-
std::size_t
get_nblocks
()¶ This member function is mostly used for ChASE with Block-Cyclic Distribution. The number of submatrix on each MPI node is
mblocks_ * nblocks_
.- Return
nblocks_
: the number of submatrices in the local matrix on each MPI node along the column direction.
-
int
get_irsrc
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
irsrc_
: the process row within 2D MPI grid over which the first row of the global matrixA
is distributed.
-
int
get_icsrc
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
icsrc_
: the process column within 2D MPI grid over which the first column of the global matrixA
is distributed.
-
std::size_t *
get_row_offs
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
r_offs_.get()
: the pointer to store offset of each subblock in local matrix along the row direction regarding the row direction of global matrixA
.
-
std::size_t *
get_row_lens
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
r_lens_.get()
: the pointer to store number of row of each subblock in local matrix along the row direction regarding the row direction of global matrixA
.
-
std::size_t *
get_row_offs_loc
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
r_offs_l_.get()
: the pointer to store offset of each subblock in local matrix along the row direction regarding the row direction of local matrix.
-
std::size_t *
get_col_offs
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
c_offs_.get()
: the pointer to store offset of each subblock of local matrix along the column direction regarding the column direction of global matrixA
.
-
std::size_t *
get_col_lens
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
c_lens_.get()
: the pointer to store length of each subblock of local matrix along the column direction regarding the column direction of global matrixA
.
-
std::size_t *
get_col_offs_loc
()¶ This member function only matters for the Block-Cyclic Distribution.
- Return
c_offs_l_.get()
: the pointer to store offset of each subblock in local matrix along the column direction regarding the column direction of local matrix.
-
void
get_offs_lens
(std::size_t *&r_offs, std::size_t *&r_lens, std::size_t *&r_offs_l, std::size_t *&c_offs, std::size_t *&c_lens, std::size_t *&c_offs_l)¶ Return the pointers to
r_offs_
,r_lens_
,r_offs_l_
,c_offs_
,c_lens_
andc_offs_l_
in single member function. This member function only matters for the Block-Cyclic Distribution.- Parameters
[in/out]
:r_offs
-> the pointer to store offset of each subblock of local matrix along the row direction regarding the global indexing of matrix A. Its size ismblocks_
, which can be obtained by get_mblocks().[in/out]
:r_lens
-> the pointer to store length of each subblock of local matrix along the row direction regarding the global indexing of matrix A. Its size ismblocks_
, which can be obtained by get_mblocks().[in/out]
:r_offs_l
-> the pointer to store offset of each subblock of local matrix along the row direction regarding the indexing of local matrix. Its size ismblocks_
, which can be obtained by get_mblocks().[in/out]
:c_offs
-> the pointer to store offset of each subblock of local matrix along the column direction regarding the global indexing of matrix A. Its size isnblocks_
, which can be obtained by get_nblocks().[in/out]
:c_lens
-> the pointer to store length of each subblock of local matrix along the column direction regarding the global indexing of matrix A. Its size isnblocks_
, which can be obtained by get_nblocks().[in/out]
:c_offs_l
-> the pointer to store offset of each subblock of local matrix along the column direction regarding the indexing of local matrix. Its size isnblocks_
, which can be obtained by get_nblocks().
-
std::string
get_dataLayout
()¶ Block-Block
Block-Cyclic
-
int *
get_coord
()¶ - Return
coord_
: the coordinates of each MPI rank in the cartesian communicator grid.
-
std::size_t
get_ldb
()¶ - Return
n_
: the leading dimension of the local matrix ofA
in the row dimension.
-
std::size_t
get_ldc
()¶ - Return
m_
: the leading dimension of the local matrix ofA
in the column dimension.
-
const std::vector<std::vector<int>> &
get_blockcounts
()¶ - Return
block_counts_
: 2D array which stores the block number of local matrix on each MPI node in each dimension.
-
const std::vector<std::vector<std::vector<int>>> &
get_blocklens
()¶ - Return
block_lens_
: 3D array which stores the size of each sublock on each MPI node in each dimension.
-
const std::vector<std::vector<std::vector<int>>> &
get_blockdispls
()¶ - Return
block_displs_
: 3D array which stores the offset of each sublock on each MPI node in each dimension.
-
const std::vector<std::vector<int>> &
get_sendlens
()¶ - Return
send_lens_
: 2D array which stores the size of local matrix on each MPI node in each dimension.
-
const std::vector<std::vector<int>> &
get_g_offsets
()¶ - Return
g_offsets_
: 2D array which stores the offset of the first block (the most up and left one) on each MPI node in each dimension regarding the global indexing of matrix A.
-
int
get_nprocs
()¶ Returns the total number of MPI nodes within MPI communicator where ChASE is working on.
- Return
nprocs_
: the total number of MPI nodes within the root MPI communicator.
-
int
get_my_rank
()¶ Returns the rank of MPI node within 2D grid.
- Return
rank_
: the rank of MPI node within 2D grid.
-
ChaseMpiMatrices<T>
create_matrices
(T *V1 = nullptr, Base<T> *ritzv = nullptr, T *V2 = nullptr, Base<T> *resid = nullptr) const¶ Create a ChaseMpiMatrices object which stores the operating matrices and vectors for ChASE-MPI.
- Parameters
V1
: am_ * max_block_
rectangular matrix.ritzv
: amax_block_
vector which stores the computed Ritz values.V2
: an_ * max_block_
rectangular matrix.resid
: amax_block_
vector which stores the residual of each computed Ritz value.
-
ChaseMpiMatrices<T>
create_matrices
(T *H, std::size_t ldh, T *V1 = nullptr, Base<T> *ritzv = nullptr, T *V2 = nullptr, Base<T> *resid = nullptr) const¶ Create a ChaseMpiMatrices object which stores the operating matrices and vectors for ChASE-MPI.
- Parameters
H
: am_ * n_
rectangular matrix, with leading dimensionldh
ldh
: leading dimension ofH
, which should be>=m_
V1
: am_ * max_block_
rectangular matrix.ritzv
: amax_block_
vector which stores the computed Ritz values.V2
: an_ * max_block_
rectangular matrix.resid
: amax_block_
vector which stores the residual of each computed Ritz value.
Private Members
-
std::size_t
N_
¶ Global size of the matrix A defining the eigenproblem.
This variable is initialized by the constructor using the value of the first of its input parameter
N
.This variable is private, it can be accessed by the member function get_N().
-
std::size_t
ldh_
¶ The leading dimension of local part of Symmetric/Hermtian matrix on each MPI proc. It should be
>= m_
.
-
std::size_t
nev_
¶ Number of desired extremal eigenpairs.
This variable is initialized by the constructor using the value of its input parameter
nev
.This variable is private, it can be accessed by the member function GetNev().
-
std::size_t
nex_
¶ Increment of the search subspace so that its total size is
nev + nex
.This variable is initialized by the constructor using the value of its input parameter
nex
.This variable is private, it can be accessed by the member function GetNex().
-
std::size_t
max_block_
¶ This variable is initialized by the constructor using the sum of its input parameter
nex
andnev
. Thus we havemax_block_=nev_+nex_
. This variable is private, it can be accessed by the member function get_max_block().
-
std::size_t
n_
¶ Column number of the local matrix on each MPI node.
For Block Distribution, this variable is initialized based on the global size of matrix
A
and dimension of the column of 2D MPI grid.For Block-Cyclic Distribution, it is determined also by block size of submatrices
nb_
, etc. This variable is private, it can be accessed by the member function get_n().
-
std::size_t
m_
¶ Row number of the local matrix on each MPI node.
For Block Distribution, this variable is initialized based on the global size of matrix
A
and dimension of the column of 2D MPI grid.For Block-Cyclic Distribution, it is determined also by block size of submatrices
mb_
, etc. This variable is private, it can be accessed by the member function get_m().
-
std::size_t
nb_
¶ Column blocking factor.
For Block Distribution, this variable equals to the variable
n_
.For Block-Cyclic Distribution, it is initialized by the parameter of the constructor
nb
, This variable is private, it can be accessed by the member function get_nb().
-
std::size_t
mb_
¶ Row blocking factor.
For Block Distribution, this variable equals to the variable
m_
.For Block-Cyclic Distribution, it is initialized by the parameter of the constructor
mb
, This variable is private, it can be accessed by the member function get_mb().
-
std::size_t
nblocks_
¶ Number of submatrices along the column direction in the local matrix. Thus each local matrix is constructed by
mblocks_ * nblocks_
blocks.For Block Distribution, this variable equals to 1.
For Block-Cyclic Distribution, it is initialized during the construction of Block-Cyclic scheme, This variable is private, it can be accessed by the member function get_nblocks().
-
std::size_t
mblocks_
¶ Number of submatrices along the row direction in the local matrix. Thus each local matrix is constructed by
mblocks_ * nblocks_
blocks.For Block Distribution, this variable equals to 1.
For Block-Cyclic Distribution, it is initialized during the construction of Block-Cyclic scheme, This variable is private, it can be accessed by the member function get_mblocks().
-
int
irsrc_
¶ Process row over which the first row of the global matrix A is distributed.
This variable matters only for the Block-Cyclic Distribution, it is initialized by the parameter irsrc of the constructor. This variable is private, it can be accessed by the member function get_irsrc().
-
int
icsrc_
¶ Process column over which the first column of the global matrix A is distributed.
This variable matters only for the Block-Cyclic Distribution, it is initialized by the parameter irsrc of the constructor. This variable is private, it can be accessed by the member function get_icsrc().
-
std::unique_ptr<std::size_t[]>
r_offs_
¶ Offset of each subblock (especially for Block-Cyclic Distribution) along the row direction regarding the global indexing of matrix
A
.This variable matters only for the Block-Cyclic Distribution, it is initialized during the setup of the Block-Cyclic Distribution scheme. This variable is private, it can be accessed by the member function get_row_offs(). The size of this array is
mblocks_
.
-
std::unique_ptr<std::size_t[]>
r_lens_
¶ Length of each subblock (especially for Block-Cyclic Distribution) along the row direction regarding the global indexing of matrix
A
.This variable matters only for the Block-Cyclic Distribution, it is initialized during the setup of the Block-Cyclic Distribution scheme. This variable is private, it can be accessed by the member function get_row_lens(). The size of this array is
mblocks_
.
-
std::unique_ptr<std::size_t[]>
r_offs_l_
¶ Offset of each subblock (especially for Block-Cyclic Distribution) along the row direction regarding the indexing of local matrix on each MPI node.
This variable matters only for the Block-Cyclic Distribution, it is initialized during the setup of the Block-Cyclic Distribution scheme. This variable is private, it can be accessed by the member function get_row_offs_loc(). The size of this array is
mblocks_
.
-
std::unique_ptr<std::size_t[]>
c_offs_
¶ Offset of each subblock (especially for Block-Cyclic Distribution) along the column direction regarding the global indexing of matrix
A
.This variable matters only for the Block-Cyclic Distribution, it is initialized during the setup of the Block-Cyclic Distribution scheme. This variable is private, it can be accessed by the member function get_col_offs(). The size of this array is
nblocks_
.
-
std::unique_ptr<std::size_t[]>
c_lens_
¶ Length of each subblock (especially for Block-Cyclic Distribution) along the column direction regarding the global indexing of matrix
A
.This variable matters only for the Block-Cyclic Distribution, it is initialized during the setup of the Block-Cyclic Distribution scheme. This variable is private, it can be accessed by the member function get_col_lens(). The size of this array is
nblocks_
.
-
std::unique_ptr<std::size_t[]>
c_offs_l_
¶ Offset of each subblock (especially for Block-Cyclic Distribution) along the column direction regarding the indexing of local matrix on each MPI node.
This variable matters only for the Block-Cyclic Distribution, it is initialized during the setup of the Block-Cyclic Distribution scheme. This variable is private, it can be accessed by the member function get_col_offs_loc(). The size of this array is
nblocks_
.
-
std::vector<std::vector<int>>
block_counts_
¶ 2D array which stores the block number of local matrix on each MPI node in each dimension.
This variable is equally shared by all the MPI nodes, which make the Block-Cyclic* scheme within each MPI node be visible to all other MPI nodes.
The variable is especically useful for the MPI communication in the case that global matrix cannot be equally distributed to each MPI node.
For example, the block number within rank
i
(this rank is within row_comm_) along the row dimension isblock_counts_[0][i]
, and the block number within rankj
(this rank is within col_comm_) along the column dimension isblock_counts_[1][j]
. For Block-Distribution, all the values in this 2D array equal to 1.
This variable is private, it can be accessed by the member function get_blockcounts().
-
std::vector<std::vector<std::vector<int>>>
block_lens_
¶ 3D array which stores the length of each sublock on each MPI node in each dimension.
This variable is equally shared by all the MPI nodes, which make the Block-Cyclic* scheme within each MPI node be visible to all other MPI nodes.
The variable is especically useful for the MPI communication in the case that global matrix cannot be equally distributed to each MPI node.
For example, the length of the block numbering
k
th within ranki
(this rank is within row_comm_) along the row dimension isblock_lens_[0][i][k]
, and the length of the block numberingk
th the block number within rankj
(this rank is within col_comm_) along the column dimension isblock_lens_[1][j][k]
.
This variable is private, it can be accessed by the member function get_blocklens().
-
std::vector<std::vector<std::vector<int>>>
block_displs_
¶ 3D array which stores the offset of each sublock on each MPI node in each dimension.
This variable is equally shared by all the MPI nodes, which make the Block-Cyclic* scheme within each MPI node be visible to all other MPI nodes.
The variable is especically useful for the MPI communication in the case that global matrix cannot be equally distributed to each MPI node.
For example, the offset of the block numbering
k
th within ranki
(this rank is within row_comm_) along the row dimension isblock_displs_[0][i][k]
, and the offset of the block numberingk
th the block number within rankj
(this rank is within col_comm_) along the column dimension isblock_displs_[1][j][k]
.
This variable is private, it can be accessed by the member function get_blockdispls().
-
std::vector<std::vector<int>>
send_lens_
¶ 2D array which stores the length of local matrix on each MPI node in each dimension.
This variable is equally shared by all the MPI nodes, which make the Block-Cyclic* scheme within each MPI node be visible to all other MPI nodes.
The variable is especically useful for the MPI communication in the case that global matrix cannot be equally distributed to each MPI node.
For example, the length of local matrix within rank
i
(this rank is within row_comm_) along the row dimension issend_lens_[0][i]
, and the block number within rankj
(this rank is within col_comm_) along the column dimension issend_lens_[1][j]
.
This variable is private, it can be accessed by the member function get_sendlens().
-
std::vector<std::vector<int>>
g_offsets_
¶ 2D array which stores the offset of the first block (the most up and left one) on each MPI node in each dimension regarding the global indexing of matrix
A
.This variable is equally shared by all the MPI nodes, which make the Block-Cyclic* scheme within each MPI node be visible to all other MPI nodes.
The variable is especically useful for the MPI communication in the case that global matrix cannot be equally distributed to each MPI node.
For example, the offset of first block of local matrix within rank
i
(this rank is within row_comm_) along the row dimension isg_offsets_[0][i]
, and the block number within rankj
(this rank is within col_comm_) along the column dimension isg_offsets_[1][j]
.
This variable is private, it can be accessed by the member function get_g_offsets().
-
MPI_Comm
comm_
¶ The MPI communicator which ChASE is working on.
This variable is initialized by the constructor using the value of its input parameters
comm
.
-
int
nprocs_
¶ Total number of MPI nodes in the MPI communicator which ChASE is working on.
-
int
rank_
¶ The rank of each MPI node within the MPI communicator which ChASE is working on.
-
std::unique_ptr<T[]>
C2_
¶ A temporary memory allocated for distributed ChASE.
This variable is initialized during the construction of ChaseMpiProperties of size
m_ * max_block_
.It can be accessed by the member function get_C2()
-
std::unique_ptr<T[]>
B2_
¶ A temporary memory allocated for distributed ChASE.
This variable is initialized during the construction of ChaseMpiProperties of size
n_ * max_block_
.It can be accessed by the member function get_B2()
-
std::unique_ptr<T[]>
A_
¶ A temporary memory allocated for distributed ChASE.
This variable is initialized during the construction of ChaseMpiProperties of size
max_block_ * max_block_
.It can be accessed by the member function get_A()
-
MPI_Comm
row_comm_
¶ The row communicator of the constructed 2D grid of MPI codes.
This variable is initialized in the constructor, after the construction of MPI 2D grid. This variable is private, it can be accessed by the member function get_row_comm().
-
MPI_Comm
col_comm_
¶ The column communicator of the constructed 2D grid of MPI codes.
This variable is initialized in the constructor, after the construction of MPI 2D grid. This variable is private, it can be accessed by the member function get_col_comm().
-
int
dims_
[2]¶ The array with two elements determines the dimension of 2D grid of MPI nodes.
For Block Distribution,
it can be initialized by
MPI_Dims_create
which creates a division of MPI ranks in a cartesian grid.it can be initialized by users with the input parameters, paramters
row_dim
andcol_dim
. More precise, we havedims_[0] = row_dim
anddims_[1] = col_dim
.
For Block-Cyclic Distribution, it is always initialized by the input paramters
row_dim
andcol_dim
. More precise, we havedims_[0] = row_dim
anddims_[1] = col_dim
.
-
int
coord_
[2]¶ The array with two elements indicates the coordinates of each MPI node within the 2D grid of MPI nodes.
This variable is determined by the properties of 2D grid of MPI node.
-
std::size_t
off_
[2]¶ The array with two elements indicates the offset of the local matrix on each MPI node regarding the global index of matrix
A
.This variable is determined by the properties of 2D grid of MPI node, the size of matrix and the scheme of distribution across the 2D grid.
-
std::string
data_layout
¶ Block-Block
Block-Cyclic