xref: /petsc/src/mat/impls/dense/mpi/mpidense.h (revision 8949adfd119cb1d44c9fb5cf3dc01f4b6a02252d)
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   Vec        lvec;                      /* local vector */
21bad5519cSBarry Smith   VecScatter Mvctx;                     /* scatter context for vector */
22ace3abfcSBarry Smith   PetscBool roworiented;                /* if true, row oriented input (default) */
23*8949adfdSHong Zhang 
24*8949adfdSHong Zhang   Mat_MatTransMatMult *atb;             /* used by MatTransposeMatMult_MPIAIJ_MPIDense */
258965ea79SLois Curfman McInnes } Mat_MPIDense;
2697304618SKris Buschelman 
275a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer);
285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
295a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
305a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*);
315a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*);
325a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
335a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
345a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
355a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetFactor_mpidense_petsc(Mat,MatFactorType,Mat*);
36