xref: /petsc/src/mat/impls/dense/mpi/mpidense.h (revision 4222ddf1d74c63827c599361b3bb0b06ad3944a0)
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*4222ddf1SHong Zhang typedef struct { /* used by MatMatMultxxx_MPIDense_MPIDense() */
7320f2790SHong Zhang   Mat            Ae,Be,Ce;           /* matrix in Elemental format */
8320f2790SHong Zhang   PetscErrorCode (*destroy)(Mat);
9320f2790SHong Zhang } Mat_MatMultDense;
10320f2790SHong Zhang 
11*4222ddf1SHong Zhang typedef struct { /* used by MatTransposeMatMultXXX_MPIDense_MPIDense() */
12c5ef1628SHong Zhang   PetscScalar    *sendbuf,*atbarray;
13d5017740SHong Zhang   PetscMPIInt    *recvcounts;
14baa3c1c6SHong Zhang   PetscErrorCode (*destroy)(Mat);
15cc48ffa7SToby Isaac   PetscMPIInt    tag;
16baa3c1c6SHong Zhang } Mat_TransMatMultDense;
17baa3c1c6SHong Zhang 
18*4222ddf1SHong Zhang typedef struct { /* used by MatMatTransposeMultxxx_MPIDense_MPIDense() */
19cc48ffa7SToby Isaac   PetscScalar    *buf[2];
20cc48ffa7SToby Isaac   PetscMPIInt    tag;
21faa55883SToby Isaac   PetscMPIInt    *recvcounts;
22faa55883SToby Isaac   PetscMPIInt    *recvdispls;
23cc48ffa7SToby Isaac   PetscErrorCode (*destroy)(Mat);
24faa55883SToby Isaac   PetscInt       alg; /* algorithm used */
25cc48ffa7SToby Isaac } Mat_MatTransMultDense;
26cc48ffa7SToby Isaac 
278965ea79SLois Curfman McInnes typedef struct {
2813f74950SBarry Smith   PetscInt    nvec;                     /* this is the n size for the vector one multiplies with */
29bad5519cSBarry Smith   Mat         A;                        /* local submatrix */
3013f74950SBarry Smith   PetscMPIInt size;                     /* size of communicator */
3113f74950SBarry Smith   PetscMPIInt rank;                     /* rank of proc in communicator */
32320f2790SHong Zhang 
338965ea79SLois Curfman McInnes   /* The following variables are used for matrix assembly */
34ace3abfcSBarry Smith   PetscBool   donotstash;               /* Flag indicationg if values should be stashed */
358965ea79SLois Curfman McInnes   MPI_Request *send_waits;              /* array of send requests */
368965ea79SLois Curfman McInnes   MPI_Request *recv_waits;              /* array of receive requests */
3713f74950SBarry Smith   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
38ea709b57SSatish Balay   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
3913f74950SBarry Smith   PetscInt    rmax;                     /* maximum message length */
408965ea79SLois Curfman McInnes 
418965ea79SLois Curfman McInnes   /* The following variables are used for matrix-vector products */
428965ea79SLois Curfman McInnes   Vec        lvec;                      /* local vector */
43bad5519cSBarry Smith   VecScatter Mvctx;                     /* scatter context for vector */
44ace3abfcSBarry Smith   PetscBool  roworiented;               /* if true, row oriented input (default) */
458949adfdSHong Zhang 
46*4222ddf1SHong Zhang   Mat_MatTransMatMult   *atb;           /* used by MatTransposeMatMultxxx_MPIAIJ_MPIDense */
47*4222ddf1SHong Zhang   Mat_TransMatMultDense *atbdense;      /* used by MatTransposeMatMultxxx_MPIDense_MPIDense */
48*4222ddf1SHong Zhang   Mat_MatMultDense      *abdense;       /* used by MatMatMultxxx_MPIDense_MPIDense */
49*4222ddf1SHong Zhang   Mat_MatTransMultDense *abtdense;      /* used by MatMatTransposeMultxxx_MPIDense_MPIDense */
508965ea79SLois Curfman McInnes } Mat_MPIDense;
5197304618SKris Buschelman 
525a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer);
535a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
547dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
555a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*);
56*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
57*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat);
58*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
595a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
60*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
61cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
62320f2790SHong Zhang 
63*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat);
64*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat);
65*4222ddf1SHong Zhang 
665242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
67320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
68*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat,Mat,PetscReal,Mat);
69*4222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat,Mat,Mat);
705242a7b1SHong Zhang #endif
71