18965ea79SLois Curfman McInnes 2c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> 38965ea79SLois Curfman McInnes 4052cd425SLois Curfman McInnes /* Data stuctures for basic parallel dense matrix */ 5052cd425SLois Curfman McInnes 6320f2790SHong Zhang typedef struct { /* used by MatMatMult_MPIDense_MPIDense() */ 7320f2790SHong Zhang Mat Ae,Be,Ce; /* matrix in Elemental format */ 8320f2790SHong Zhang PetscErrorCode (*destroy)(Mat); 9320f2790SHong Zhang } Mat_MatMultDense; 10320f2790SHong Zhang 11baa3c1c6SHong Zhang typedef struct { /* used by MatTransposeMatMult_MPIDense_MPIDense() */ 12c5ef1628SHong Zhang PetscScalar *sendbuf,*atbarray; 13d5017740SHong Zhang PetscMPIInt *recvcounts; 14baa3c1c6SHong Zhang PetscErrorCode (*destroy)(Mat); 15*cc48ffa7SToby Isaac PetscMPIInt tag; 16baa3c1c6SHong Zhang } Mat_TransMatMultDense; 17baa3c1c6SHong Zhang 18*cc48ffa7SToby Isaac typedef struct { /* used by MatMatTransposeMult_MPIDense_MPIDense() */ 19*cc48ffa7SToby Isaac PetscScalar *buf[2]; 20*cc48ffa7SToby Isaac PetscMPIInt tag; 21*cc48ffa7SToby Isaac PetscErrorCode (*destroy)(Mat); 22*cc48ffa7SToby Isaac } Mat_MatTransMultDense; 23*cc48ffa7SToby Isaac 248965ea79SLois Curfman McInnes typedef struct { 2513f74950SBarry Smith PetscInt nvec; /* this is the n size for the vector one multiplies with */ 26bad5519cSBarry Smith Mat A; /* local submatrix */ 2713f74950SBarry Smith PetscMPIInt size; /* size of communicator */ 2813f74950SBarry Smith PetscMPIInt rank; /* rank of proc in communicator */ 29320f2790SHong Zhang 308965ea79SLois Curfman McInnes /* The following variables are used for matrix assembly */ 31ace3abfcSBarry Smith PetscBool donotstash; /* Flag indicationg if values should be stashed */ 328965ea79SLois Curfman McInnes MPI_Request *send_waits; /* array of send requests */ 338965ea79SLois Curfman McInnes MPI_Request *recv_waits; /* array of receive requests */ 3413f74950SBarry Smith PetscInt nsends,nrecvs; /* numbers of sends and receives */ 35ea709b57SSatish Balay PetscScalar *svalues,*rvalues; /* sending and receiving data */ 3613f74950SBarry Smith PetscInt rmax; /* maximum message length */ 378965ea79SLois Curfman McInnes 388965ea79SLois Curfman McInnes /* The following variables are used for matrix-vector products */ 398965ea79SLois Curfman McInnes Vec lvec; /* local vector */ 40bad5519cSBarry Smith VecScatter Mvctx; /* scatter context for vector */ 41ace3abfcSBarry Smith PetscBool roworiented; /* if true, row oriented input (default) */ 428949adfdSHong Zhang 438949adfdSHong Zhang Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult_MPIAIJ_MPIDense */ 44baa3c1c6SHong Zhang Mat_TransMatMultDense *atbdense; /* used by MatTransposeMatMult_MPIDense_MPIDense */ 45320f2790SHong Zhang Mat_MatMultDense *abdense; /* used by MatMatMult_MPIDense_MPIDense */ 46*cc48ffa7SToby Isaac Mat_MatTransMultDense *abtdense; /* used by MatMatTransposeMult_MPIDense_MPIDense */ 478965ea79SLois Curfman McInnes } Mat_MPIDense; 4897304618SKris Buschelman 495a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 505a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 517dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 525a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*); 535a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 545a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 555a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 565a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 57cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 58cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 59cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 60320f2790SHong Zhang 615242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL) 62320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 63320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 645242a7b1SHong Zhang #endif 65