Go to the documentation of this file. 1 #ifndef PARGEMSLR_PARALLEL_CSR_MATRIX_H
2 #define PARGEMSLR_PARALLEL_CSR_MATRIX_H
10 #include "../utils/utils.hpp"
11 #include "../utils/structs.hpp"
52 std::vector<void*> _send_buff_v2;
58 std::vector<void*> _recv_buff_v2;
66 int _buff_d_unit_size;
72 std::vector<void*> _send_buff_v2_d;
78 std::vector<void*> _recv_buff_v2_d;
208 #ifdef PARGEMSLR_CUDA
216 int CreateDeviceBuffer(
int unitsize);
250 template <
typename T>
261 template <
typename T>
271 template <
typename T>
283 template <
typename T,
class VectorType>
284 int DataTransfer(
const std::vector<VectorType> &vec_in, std::vector<VectorType> &vec_out,
int loc_in,
int loc_out);
294 template <
typename T,
class VectorType>
304 template <
typename T,
class VectorType>
316 template <
typename T>
327 template <
typename T>
337 template <
typename T>
346 template <
typename T>
357 long int _nrow_global;
363 long int _ncol_global;
369 long int _nrow_start;
375 long int _ncol_start;
394 bool _is_offd_map_sorted;
446 typename std::conditional<PargemslrIsDoublePrecision<T>::value,
529 int Setup(
int nrow_local,
long int nrow_start,
long int nrow_global,
int ncol_local,
long int ncol_start,
long int ncol_global,
parallel_log &parlog);
667 virtual int Fill(
const T &v);
675 virtual int MoveData(
const int &location);
683 virtual int Scale(
const T &alpha);
815 int Laplacian(
int nx,
int ny,
int nz,
int dx,
int dy,
int dz, T alphax, T alphay, T alphaz, T shift,
parallel_log &parlog,
bool rand_perturb =
false);
852 int LaplacianWithPartition(
int nx,
int ny,
int nz,
int dx,
int dy,
int dz,
int d2x,
int d2y,
int d2z, T alphax, T alphay, T alphaz, T shift,
parallel_log &parlog,
bool rand_perturb =
false);
895 template <
typename T1>
899 return PARGEMSLR_SUCCESS;
908 template <
typename T1>
912 return PARGEMSLR_SUCCESS;
922 return this->_separator_ndom;
932 return this->_separator_domi;
937 typedef ParallelCsrMatrixClass<float> matrix_csr_par_float;
938 typedef ParallelCsrMatrixClass<double> matrix_csr_par_double;
939 typedef ParallelCsrMatrixClass<complexs> matrix_csr_par_complexs;
940 typedef ParallelCsrMatrixClass<complexd> matrix_csr_par_complexd;
int SetComplexShift(T1 diagonal_shift)
Set the diagonal complex shift for some preconditioner options.
Definition: parallel_csr_matrix.hpp:896
int SetOffdMatSorted(bool is_offd_map_sorted)
Set the off-diagonal map vector to sorted/unsorted.
Definition: parallel_csr_matrix.cpp:1529
int GetDiagScale(T &scale)
Get the value leads to max diagonal no.
Definition: parallel_csr_matrix.cpp:1710
ParallelCsrMatrixClass()
The constructor of ParallelCsrMatrixClass.
Definition: parallel_csr_matrix.cpp:1115
virtual long int GetNumNonzeros() const
Get the number of nonzeros in this matrix.
Definition: parallel_csr_matrix.cpp:1416
bool _is_waiting
Is true when there is communication on going.
Definition: parallel_csr_matrix.hpp:143
int Helmholtz(int n, int dx, int dy, int dz, T w, parallel_log &parlog)
Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = f.
Definition: parallel_csr_matrix.cpp:3735
int GetComplexShift(T1 &diagonal_shift)
Get the diagonal complex shift for some preconditioner options.
Definition: csr_matrix.cpp:1723
int SetupMatvec()
Setup the matvec_helper for parallel matvec.
Definition: parallel_csr_matrix.cpp:1731
virtual int Scale(const T &alpha)
Scale the current csr matrix.
Definition: parallel_csr_matrix.cpp:1696
virtual int MatVec(char trans, const T &alpha, const VectorClass< T > &x, const T &beta, VectorClass< T > &y)
In place parallel csr Matrix-Vector product ==> y := alpha*A*x + beta*y, or y := alpha*A'*x + beta*y.
Definition: parallel_csr_matrix.cpp:2105
int DataTransferStart(const VectorVirtualClass< T > &vec_in, int loc_in)
Start communication using the buffer inside this comm helper.
Definition: parallel_csr_matrix.cpp:447
int Setup(int nrow_local, int ncol_local, parallel_log &parlog)
Free the current parallel matrix, and set up new paralle matrix.
Definition: parallel_csr_matrix.cpp:1295
int LaplacianWithPartition(int nx, int ny, int nz, int dx, int dy, int dz, int d2x, int d2y, int d2z, T alphax, T alphay, T alphaz, T shift, parallel_log &parlog, bool rand_perturb=false)
Generate Laplacian matrix on the host memory, 5-pt for 2D problem and 7-pt for 3D problem.
Definition: parallel_csr_matrix.cpp:3068
CSR matrix data structure.
int MoveData(int location)
Move the index data to other location in the memory.
Definition: parallel_csr_matrix.cpp:254
Tell if a value is a parallel data structure.
Definition: utils.hpp:444
int MoveRecvData(int location)
Move the recv data to other location in the memory.
Definition: parallel_csr_matrix.cpp:284
vector< MPI_Request > _requests_v
Vector of MPI_Request for communication.
Definition: parallel_csr_matrix.hpp:131
int GetComplexShift(T1 &diagonal_shift)
Get the diagonal complex shift for some preconditioner options.
Definition: parallel_csr_matrix.hpp:909
virtual int GetNumColsLocal() const
Get the local number of columns of the matrix.
Definition: parallel_csr_matrix.cpp:1366
int SetComplexShift(T1 diagonal_shift)
Set the diagonal complex shift for some preconditioner options.
Definition: csr_matrix.hpp:666
long int GetColStartGlobal() const
Get the global number of the first column of the matrix.
Definition: parallel_csr_matrix.cpp:1406
vector_int _recv_from_v
The vector stores the MPI ranks that myid needs to recv external vector from.
Definition: parallel_csr_matrix.hpp:106
vector_long & GetOffdMap()
Get the off-diagonal column map vector.
Definition: parallel_csr_matrix.cpp:1509
CsrMatrixClass< T > & GetOffdMat()
Get the off-diagonal matrix.
Definition: parallel_csr_matrix.cpp:1499
int MatMat(const T &alpha, const ParallelCsrMatrixClass< T > &A, char transa, const ParallelCsrMatrixClass< T > &B, char transb, const T &beta)
Dense Matrix-Matrix multiplication, C = alpha*transa(A)*transb(B)+beta*C where C is this matrix.
Definition: parallel_csr_matrix.cpp:2233
The virtual class of real/complex vector class.
Definition: vector.hpp:126
ParallelCsrMatrixClass< T > & operator=(const ParallelCsrMatrixClass< T > &mat)
The = operator of ParallelCsrMatrixClass.
Definition: parallel_csr_matrix.cpp:1190
The class of parallel real/complex vector.
Definition: structs.hpp:14
virtual ~ParallelCsrMatrixClass()
The destructor of ParallelCsrMatrixClass.
Definition: parallel_csr_matrix.cpp:1254
CommunicationHelperClass & operator=(const CommunicationHelperClass &comm_helper)
The = operator of CommunicationHelperClass.
Definition: parallel_csr_matrix.cpp:104
The class of sequential real/complex vector.
Definition: structs.hpp:12
bool _is_nbhd_built
If the neiborhood information (send_to and recv_from) is set in advance.
Definition: parallel_csr_matrix.hpp:149
long int GetNumColsGlobal() const
Get the global number of columns of the matrix.
Definition: parallel_csr_matrix.cpp:1386
The data structure for parallel matvec helper, store communication inforamtion for parallel matvec,...
Definition: parallel_csr_matrix.hpp:27
int SetupMatvecOver()
Finishing up the setup phase of matvec_helper for parallel matvec.
Definition: parallel_csr_matrix.cpp:2060
int DataTransferStartReverse(const VectorVirtualClass< T > &vec_in, int loc_in)
Start communication using the buffer inside this comm helper in the opposite direction.
Definition: parallel_csr_matrix.cpp:853
int GetGraphArrays(vector_long &vtxdist, vector_long &xadj, vector_long &adjncy)
Get the vector arrays for the graph, used as inputs for ParMETIS. Diagonal entry is removed.
Definition: parallel_csr_matrix.cpp:1435
Class of parallel CSR matrices.
Definition: structs.hpp:18
The virtual matrix classes.
virtual int GetNumRowsLocal() const
Get the local number of rows of the matrix.
Definition: parallel_csr_matrix.cpp:1356
vector_int _send_to_v
The vector stores the MPI ranks that myid needs to send local vector to.
Definition: parallel_csr_matrix.hpp:100
long int GetNumRowsGlobal() const
Get the global number of rows of the matrix.
Definition: parallel_csr_matrix.cpp:1376
int Laplacian(int nx, int ny, int nz, int dx, int dy, int dz, T alphax, T alphay, T alphaz, T shift, parallel_log &parlog, bool rand_perturb=false)
Generate Laplacian matrix on the host memory, 5-pt for 2D problem and 7-pt for 3D problem.
Definition: parallel_csr_matrix.cpp:2519
The template class complex.
Definition: complex.hpp:24
int SubMatrix(vector_long &rows, vector_long &cols, int location, ParallelCsrMatrixClass< T > &parcsrmat_out)
Extract parallel submatrix from parallel CSR matrix.
Definition: parallel_csr_matrix.cpp:2244
CommunicationHelperClass _comm_helper
The data struxture for the parallel matvec. Stores communication inforamtion, and working buffers.
Definition: parallel_csr_matrix.hpp:454
int & GetSeparatorNumSubdomains()
Get number of subdomains on the top level if setted in advange.
Definition: parallel_csr_matrix.hpp:920
long int GetRowStartGlobal() const
Get the global number of the first row of the matrix.
Definition: parallel_csr_matrix.cpp:1396
int HelmholtzWithPartition(int n, int dx, int dy, int dz, int d2x, int d2y, int d2z, T w, parallel_log &parlog)
Generate 3D Helmholtz matrix with 7-pt on [0,1]^3. -\Delta u - w^2 u = f.
Definition: parallel_csr_matrix.cpp:4335
int PlotPatternGnuPlot(const char *datafilename, int pttype=0)
Plot the pattern of the parallel csr matrix using gnuplot. Similar to spy in the MATLAB....
Definition: parallel_csr_matrix.cpp:5384
virtual int Eye()
Create an indentity matrix.
Definition: parallel_csr_matrix.cpp:1638
The virtual class of all matrix classes.
Definition: matrix.hpp:23
int _n_in
The size of the input vector.
Definition: parallel_csr_matrix.hpp:88
class pargemslr::ParallelLogClass parallel_log
The data structure for parallel computing, including data structures for MPI and CUDA.
virtual int Fill(const T &v)
Fill the matrix with constant value.
Definition: parallel_csr_matrix.cpp:1683
bool _is_ready
Is true when the comm_helper ready to use.
Definition: parallel_csr_matrix.hpp:137
int Transpose()
Transpost the current matrix.
Definition: parallel_csr_matrix.cpp:2479
CommunicationHelperClass _trans_comm_helper
The data struxture for the parallel matvec of A transpose. Stores communication inforamtion,...
Definition: parallel_csr_matrix.hpp:461
int DataTransferOverReverse(VectorVirtualClass< T > &vec_out, int loc_out)
Finish communication using the buffer inside this comm helper in the opposite direction.
Definition: parallel_csr_matrix.cpp:997
std::vector< vector_int > _send_idx_v2
2D vector, send_idx_v2[i] is the index in local vector that need to be sent to proc....
Definition: parallel_csr_matrix.hpp:112
std::vector< vector_int > _idx_helper_v2
Helper array for building the comm helper.
Definition: parallel_csr_matrix.hpp:125
int SetupMatvecStart()
Setup the matvec_helper for parallel matvec.
Definition: parallel_csr_matrix.cpp:1746
virtual int MoveData(const int &location)
Move the data to another memory location.
Definition: parallel_csr_matrix.cpp:1658
std::conditional< PargemslrIsDoublePrecision< T >::value, double, float >::type _diagonal_shift
The diagonal complex shift in some preconditioners. The preconditioner would be built on A+_diagonal_...
Definition: parallel_csr_matrix.hpp:448
CsrMatrixClass< T > & GetDiagMat()
Get the diagonal matrix.
Definition: parallel_csr_matrix.cpp:1489
int SortOffdMap()
Sort the off-diagonal map vector.
Definition: parallel_csr_matrix.cpp:1540
virtual int GetDataLocation() const
Get the data location of the matrix.
Definition: parallel_csr_matrix.cpp:1345
int SetupVectorPtrStr(ParallelVectorClass< T > &vec)
Update the structure of a parallel vector to have same distribution as the row distribution of this m...
Definition: parallel_csr_matrix.cpp:1628
int DataTransfer(const VectorVirtualClass< T > &vec_in, VectorVirtualClass< T > &vec_out, int loc_in, int loc_out)
Apply communication using the buffer inside this comm helper.
Definition: parallel_csr_matrix.cpp:433
int _n_out
The size of the output vector.
Definition: parallel_csr_matrix.hpp:94
The virtual class of all vector classes.
Definition: vector.hpp:22
Tell if a value is a complex value.
Definition: complex.hpp:684
int MatVecOffd(char trans, const T &alpha, const VectorClass< T > &x, const T &beta, VectorClass< T > &y)
In place parallel csr Matrix-Vector product with offdiagonal block only ==> y := alpha*A*x + beta*y,...
Definition: parallel_csr_matrix.cpp:2195
int MoveSendData(int location)
Move the send data to other location in the memory.
Definition: parallel_csr_matrix.cpp:264
CommunicationHelperClass()
The constructor of CommunicationHelperClass.
Definition: parallel_csr_matrix.cpp:25
int CreateHostBuffer(int unitsize)
Create buffer.
Definition: parallel_csr_matrix.cpp:304
int Clear()
Free the current matvec helper, set everything to 0.
Definition: parallel_csr_matrix.cpp:176
~CommunicationHelperClass()
The destructor of CommunicationHelperClass.
Definition: parallel_csr_matrix.cpp:171
int DataTransferReverse(const VectorVirtualClass< T > &vec_in, VectorVirtualClass< T > &vec_out, int loc_in, int loc_out)
Apply communication using the buffer inside this comm helper in the opposite direction.
Definition: parallel_csr_matrix.cpp:839
Class of CSR/CSC matrices.
Definition: structs.hpp:16
std::vector< vector_int > _recv_idx_v2
2D vector, recv_idx_v2[i] is the index in external vector that need to be received from proc....
Definition: parallel_csr_matrix.hpp:119
bool & IsOffdMapSorted()
Tell if the off-diagonal map vector is sorted.
Definition: parallel_csr_matrix.cpp:1519
vector_int & GetSeparatorDomi()
Number of subdomain for each node, -1 is the separator.
Definition: parallel_csr_matrix.hpp:930
virtual int Clear()
Free the current matrix.
Definition: parallel_csr_matrix.cpp:1264
int DataTransferOver(VectorVirtualClass< T > &vec_out, int loc_out)
Finish communication using the buffer inside this comm helper.
Definition: parallel_csr_matrix.cpp:591
int ReadFromSingleMMFile(const char *matfile, int idxin, parallel_log &parlog)
Read matrix on the host memory, matrix market format.
Definition: parallel_csr_matrix.cpp:5135