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*baa3c1c6SHong Zhang typedef struct { /* used by MatTransposeMatMult_MPIDense_MPIDense() */ 7*baa3c1c6SHong Zhang PetscScalar *sendbuf,*recvbuf; 8*baa3c1c6SHong Zhang PetscInt *recvcounts; 9*baa3c1c6SHong Zhang Mat C_loc; 10*baa3c1c6SHong Zhang PetscErrorCode (*destroy)(Mat); 11*baa3c1c6SHong Zhang } Mat_TransMatMultDense; 12*baa3c1c6SHong Zhang 138965ea79SLois Curfman McInnes typedef struct { 1413f74950SBarry Smith PetscInt nvec; /* this is the n size for the vector one multiplies with */ 15bad5519cSBarry Smith Mat A; /* local submatrix */ 1613f74950SBarry Smith PetscMPIInt size; /* size of communicator */ 1713f74950SBarry Smith PetscMPIInt rank; /* rank of proc in communicator */ 188965ea79SLois Curfman McInnes /* The following variables are used for matrix assembly */ 19ace3abfcSBarry Smith PetscBool donotstash; /* Flag indicationg if values should be stashed */ 208965ea79SLois Curfman McInnes MPI_Request *send_waits; /* array of send requests */ 218965ea79SLois Curfman McInnes MPI_Request *recv_waits; /* array of receive requests */ 2213f74950SBarry Smith PetscInt nsends,nrecvs; /* numbers of sends and receives */ 23ea709b57SSatish Balay PetscScalar *svalues,*rvalues; /* sending and receiving data */ 2413f74950SBarry Smith PetscInt rmax; /* maximum message length */ 258965ea79SLois Curfman McInnes 268965ea79SLois Curfman McInnes /* The following variables are used for matrix-vector products */ 278965ea79SLois Curfman McInnes Vec lvec; /* local vector */ 28bad5519cSBarry Smith VecScatter Mvctx; /* scatter context for vector */ 29ace3abfcSBarry Smith PetscBool roworiented; /* if true, row oriented input (default) */ 308949adfdSHong Zhang 318949adfdSHong Zhang Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult_MPIAIJ_MPIDense */ 32*baa3c1c6SHong Zhang Mat_TransMatMultDense *atbdense; /* used by MatTransposeMatMult_MPIDense_MPIDense */ 338965ea79SLois Curfman McInnes } Mat_MPIDense; 3497304618SKris Buschelman 355a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 365a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 375a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 385a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*); 395a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 405a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 415a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 425a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 43cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 44cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 45cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 46