xref: /petsc/src/mat/impls/dense/mpi/mpidense.h (revision 6947451f174971b5cbba4c8ab49dd392302c23da)
18965ea79SLois Curfman McInnes 
2c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h>
3637a0070SStefano Zampini #include <petscsf.h>
48965ea79SLois Curfman McInnes 
5052cd425SLois Curfman McInnes /*  Data stuctures for basic parallel dense matrix  */
6052cd425SLois Curfman McInnes 
74222ddf1SHong Zhang typedef struct { /* used by MatMatMultxxx_MPIDense_MPIDense() */
8320f2790SHong Zhang   Mat            Ae,Be,Ce;           /* matrix in Elemental format */
9320f2790SHong Zhang   PetscErrorCode (*destroy)(Mat);
10320f2790SHong Zhang } Mat_MatMultDense;
11320f2790SHong Zhang 
124222ddf1SHong Zhang typedef struct { /* used by MatTransposeMatMultXXX_MPIDense_MPIDense() */
13637a0070SStefano Zampini   PetscScalar    *sendbuf;
14637a0070SStefano Zampini   Mat            atb;
15d5017740SHong Zhang   PetscMPIInt    *recvcounts;
16baa3c1c6SHong Zhang   PetscErrorCode (*destroy)(Mat);
17cc48ffa7SToby Isaac   PetscMPIInt    tag;
18baa3c1c6SHong Zhang } Mat_TransMatMultDense;
19baa3c1c6SHong Zhang 
204222ddf1SHong Zhang typedef struct { /* used by MatMatTransposeMultxxx_MPIDense_MPIDense() */
21cc48ffa7SToby Isaac   PetscScalar    *buf[2];
22cc48ffa7SToby Isaac   PetscMPIInt    tag;
23faa55883SToby Isaac   PetscMPIInt    *recvcounts;
24faa55883SToby Isaac   PetscMPIInt    *recvdispls;
25cc48ffa7SToby Isaac   PetscErrorCode (*destroy)(Mat);
26faa55883SToby Isaac   PetscInt       alg; /* algorithm used */
27cc48ffa7SToby Isaac } Mat_MatTransMultDense;
28cc48ffa7SToby Isaac 
298965ea79SLois Curfman McInnes typedef struct {
30bad5519cSBarry Smith   Mat         A;                        /* local submatrix */
3113f74950SBarry Smith   PetscMPIInt size;                     /* size of communicator */
3213f74950SBarry Smith   PetscMPIInt rank;                     /* rank of proc in communicator */
33320f2790SHong Zhang 
348965ea79SLois Curfman McInnes   /* The following variables are used for matrix assembly */
35ace3abfcSBarry Smith   PetscBool   donotstash;               /* Flag indicationg if values should be stashed */
368965ea79SLois Curfman McInnes   MPI_Request *send_waits;              /* array of send requests */
378965ea79SLois Curfman McInnes   MPI_Request *recv_waits;              /* array of receive requests */
3813f74950SBarry Smith   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
39ea709b57SSatish Balay   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
4013f74950SBarry Smith   PetscInt    rmax;                     /* maximum message length */
418965ea79SLois Curfman McInnes 
428965ea79SLois Curfman McInnes   /* The following variables are used for matrix-vector products */
438965ea79SLois Curfman McInnes   Vec        lvec;                      /* local vector */
44637a0070SStefano Zampini   PetscSF    Mvctx;                     /* for mat-mult communications */
45ace3abfcSBarry Smith   PetscBool  roworiented;               /* if true, row oriented input (default) */
468949adfdSHong Zhang 
47*6947451fSStefano Zampini   /* Support for MatDenseGetColumnVec */
48*6947451fSStefano Zampini   Vec               cvec;      /* vector representation of a given column */
49*6947451fSStefano Zampini   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
50*6947451fSStefano Zampini   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
51*6947451fSStefano Zampini 
524222ddf1SHong Zhang   Mat_MatTransMatMult   *atb;           /* used by MatTransposeMatMultxxx_MPIAIJ_MPIDense */
534222ddf1SHong Zhang   Mat_TransMatMultDense *atbdense;      /* used by MatTransposeMatMultxxx_MPIDense_MPIDense */
544222ddf1SHong Zhang   Mat_MatMultDense      *abdense;       /* used by MatMatMultxxx_MPIDense_MPIDense */
554222ddf1SHong Zhang   Mat_MatTransMultDense *abtdense;      /* used by MatMatTransposeMultxxx_MPIDense_MPIDense */
568965ea79SLois Curfman McInnes } Mat_MPIDense;
5797304618SKris Buschelman 
585a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer);
595a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
607dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
615a576424SJed Brown PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*);
624222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
634222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat);
644222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
655a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
664222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
67cb20be35SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
68320f2790SHong Zhang 
694222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat);
704222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat);
714222ddf1SHong Zhang 
725242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
73320f2790SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
744222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat,Mat,PetscReal,Mat);
754222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat,Mat,Mat);
765242a7b1SHong Zhang #endif
77