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 6*320f2790SHong Zhang typedef struct { /* used by MatMatMult_MPIDense_MPIDense() */ 7*320f2790SHong Zhang Mat Ae,Be,Ce; /* matrix in Elemental format */ 8*320f2790SHong Zhang PetscErrorCode (*destroy)(Mat); 9*320f2790SHong Zhang } Mat_MatMultDense; 10*320f2790SHong Zhang 11baa3c1c6SHong Zhang typedef struct { /* used by MatTransposeMatMult_MPIDense_MPIDense() */ 12c5ef1628SHong Zhang PetscScalar *sendbuf,*atbarray; 13d5017740SHong Zhang PetscMPIInt *recvcounts; 14baa3c1c6SHong Zhang PetscErrorCode (*destroy)(Mat); 15baa3c1c6SHong Zhang } Mat_TransMatMultDense; 16baa3c1c6SHong Zhang 178965ea79SLois Curfman McInnes typedef struct { 1813f74950SBarry Smith PetscInt nvec; /* this is the n size for the vector one multiplies with */ 19bad5519cSBarry Smith Mat A; /* local submatrix */ 2013f74950SBarry Smith PetscMPIInt size; /* size of communicator */ 2113f74950SBarry Smith PetscMPIInt rank; /* rank of proc in communicator */ 22*320f2790SHong Zhang 238965ea79SLois Curfman McInnes /* The following variables are used for matrix assembly */ 24ace3abfcSBarry Smith PetscBool donotstash; /* Flag indicationg if values should be stashed */ 258965ea79SLois Curfman McInnes MPI_Request *send_waits; /* array of send requests */ 268965ea79SLois Curfman McInnes MPI_Request *recv_waits; /* array of receive requests */ 2713f74950SBarry Smith PetscInt nsends,nrecvs; /* numbers of sends and receives */ 28ea709b57SSatish Balay PetscScalar *svalues,*rvalues; /* sending and receiving data */ 2913f74950SBarry Smith PetscInt rmax; /* maximum message length */ 308965ea79SLois Curfman McInnes 318965ea79SLois Curfman McInnes /* The following variables are used for matrix-vector products */ 328965ea79SLois Curfman McInnes Vec lvec; /* local vector */ 33bad5519cSBarry Smith VecScatter Mvctx; /* scatter context for vector */ 34ace3abfcSBarry Smith PetscBool roworiented; /* if true, row oriented input (default) */ 358949adfdSHong Zhang 368949adfdSHong Zhang Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult_MPIAIJ_MPIDense */ 37baa3c1c6SHong Zhang Mat_TransMatMultDense *atbdense; /* used by MatTransposeMatMult_MPIDense_MPIDense */ 38*320f2790SHong Zhang Mat_MatMultDense *abdense; /* used by MatMatMult_MPIDense_MPIDense */ 398965ea79SLois Curfman McInnes } Mat_MPIDense; 4097304618SKris Buschelman 415a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 425a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 435a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 445a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*); 455a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 465a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 475a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 485a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 49cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 50cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 51cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 52*320f2790SHong Zhang 53*320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 54*320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 55*320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 56