18965ea79SLois Curfman McInnes 2c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> 3637a0070SStefano Zampini #include <petscsf.h> 48965ea79SLois Curfman McInnes 5da81f932SPierre Jolivet /* Data structures for basic parallel dense matrix */ 6052cd425SLois Curfman McInnes 74222ddf1SHong Zhang typedef struct { /* used by MatMatMultxxx_MPIDense_MPIDense() */ 8320f2790SHong Zhang Mat Ae, Be, Ce; /* matrix in Elemental format */ 9320f2790SHong Zhang } Mat_MatMultDense; 10320f2790SHong Zhang 114222ddf1SHong Zhang typedef struct { /* used by MatTransposeMatMultXXX_MPIDense_MPIDense() */ 12637a0070SStefano Zampini PetscScalar *sendbuf; 13637a0070SStefano Zampini Mat atb; 14d5017740SHong Zhang PetscMPIInt *recvcounts; 15cc48ffa7SToby Isaac PetscMPIInt tag; 16baa3c1c6SHong Zhang } Mat_TransMatMultDense; 17baa3c1c6SHong Zhang 184222ddf1SHong Zhang typedef struct { /* used by MatMatTransposeMultxxx_MPIDense_MPIDense() */ 19cc48ffa7SToby Isaac PetscScalar *buf[2]; 20cc48ffa7SToby Isaac PetscMPIInt tag; 21faa55883SToby Isaac PetscMPIInt *recvcounts; 22faa55883SToby Isaac PetscMPIInt *recvdispls; 23faa55883SToby Isaac PetscInt alg; /* algorithm used */ 24cc48ffa7SToby Isaac } Mat_MatTransMultDense; 25cc48ffa7SToby Isaac 268965ea79SLois Curfman McInnes typedef struct { 27bad5519cSBarry Smith Mat A; /* local submatrix */ 28320f2790SHong Zhang 298965ea79SLois Curfman McInnes /* The following variables are used for matrix assembly */ 30ace3abfcSBarry Smith PetscBool donotstash; /* Flag indicationg if values should be stashed */ 318965ea79SLois Curfman McInnes MPI_Request *send_waits; /* array of send requests */ 328965ea79SLois Curfman McInnes MPI_Request *recv_waits; /* array of receive requests */ 3313f74950SBarry Smith PetscInt nsends, nrecvs; /* numbers of sends and receives */ 34ea709b57SSatish Balay PetscScalar *svalues, *rvalues; /* sending and receiving data */ 3513f74950SBarry Smith PetscInt rmax; /* maximum message length */ 368965ea79SLois Curfman McInnes 378965ea79SLois Curfman McInnes /* The following variables are used for matrix-vector products */ 388965ea79SLois Curfman McInnes Vec lvec; /* local vector */ 39637a0070SStefano Zampini PetscSF Mvctx; /* for mat-mult communications */ 40ace3abfcSBarry Smith PetscBool roworiented; /* if true, row oriented input (default) */ 418949adfdSHong Zhang 425ea7661aSPierre Jolivet /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */ 435ea7661aSPierre Jolivet Mat cmat; /* matrix representation of a given subset of columns */ 446947451fSStefano Zampini Vec cvec; /* vector representation of a given column */ 456947451fSStefano Zampini const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */ 466947451fSStefano Zampini PetscInt vecinuse; /* if cvec is in use (col = vecinuse-1) */ 475ea7661aSPierre Jolivet PetscInt matinuse; /* if cmat is in use (cbegin = matinuse-1) */ 488965ea79SLois Curfman McInnes } Mat_MPIDense; 4997304618SKris Buschelman 505a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 517dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]); 524222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat); 53320f2790SHong Zhang 544222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat); 554222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat); 564222ddf1SHong Zhang 575242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 584222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat, Mat, PetscReal, Mat); 594222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat, Mat, Mat); 605242a7b1SHong Zhang #endif 614e29119aSPierre Jolivet PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype, MatReuse reuse, Mat *newmat); 62*4742e46bSJacob Faibussowitsch 63*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatShift_MPIDense(Mat, PetscScalar); 64*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_MPIDense(Mat, PetscInt, Vec *); 65*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_MPIDense(Mat, PetscInt, Vec *); 66*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_MPIDense(Mat, PetscInt, Vec *); 67*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_MPIDense(Mat, PetscInt, Vec *); 68*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_MPIDense(Mat, PetscInt, Vec *); 69*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_MPIDense(Mat, PetscInt, Vec *); 70*4742e46bSJacob Faibussowitsch 71*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatCreate_MPIDense(Mat); 72*4742e46bSJacob Faibussowitsch 73*4742e46bSJacob Faibussowitsch #if PetscDefined(HAVE_CUDA) 74*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_MPIDense_MPIDenseCUDA(Mat, MatType, MatReuse, Mat *); 75*4742e46bSJacob Faibussowitsch #endif 76*4742e46bSJacob Faibussowitsch 77*4742e46bSJacob Faibussowitsch #if PetscDefined(HAVE_HIP) 78*4742e46bSJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatConvert_MPIDense_MPIDenseHIP(Mat, MatType, MatReuse, Mat *); 79*4742e46bSJacob Faibussowitsch #endif 80