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 68965ea79SLois Curfman McInnes typedef struct { 713f74950SBarry Smith PetscInt nvec; /* this is the n size for the vector one multiplies with */ 8bad5519cSBarry Smith Mat A; /* local submatrix */ 913f74950SBarry Smith PetscMPIInt size; /* size of communicator */ 1013f74950SBarry Smith PetscMPIInt rank; /* rank of proc in communicator */ 118965ea79SLois Curfman McInnes /* The following variables are used for matrix assembly */ 12ace3abfcSBarry Smith PetscBool donotstash; /* Flag indicationg if values should be stashed */ 138965ea79SLois Curfman McInnes MPI_Request *send_waits; /* array of send requests */ 148965ea79SLois Curfman McInnes MPI_Request *recv_waits; /* array of receive requests */ 1513f74950SBarry Smith PetscInt nsends,nrecvs; /* numbers of sends and receives */ 16ea709b57SSatish Balay PetscScalar *svalues,*rvalues; /* sending and receiving data */ 1713f74950SBarry Smith PetscInt rmax; /* maximum message length */ 188965ea79SLois Curfman McInnes 198965ea79SLois Curfman McInnes /* The following variables are used for matrix-vector products */ 208965ea79SLois Curfman McInnes 218965ea79SLois Curfman McInnes Vec lvec; /* local vector */ 22bad5519cSBarry Smith VecScatter Mvctx; /* scatter context for vector */ 23052cd425SLois Curfman McInnes 24ace3abfcSBarry Smith PetscBool roworiented; /* if true, row oriented input (default) */ 258965ea79SLois Curfman McInnes } Mat_MPIDense; 2697304618SKris Buschelman 27*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 28*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 29*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 30*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*); 31*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 32*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 33*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 34*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 35*5a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetFactor_mpidense_petsc(Mat,MatFactorType,Mat*); 3601b82886SBarry Smith 3701b82886SBarry Smith 38