xref: /petsc/src/mat/impls/dense/mpi/mpidense.h (revision 5ea7661a5c6ebfe94618efd9ae12232cbd7a6c1a)
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 } Mat_MatMultDense;
10320f2790SHong Zhang 
114222ddf1SHong Zhang typedef struct { /* used by MatTransposeMatMultXXX_MPIDense_MPIDense() */
12637a0070SStefano Zampini   PetscScalar *sendbuf;
13637a0070SStefano Zampini   Mat         atb;
14d5017740SHong Zhang   PetscMPIInt *recvcounts;
15cc48ffa7SToby Isaac   PetscMPIInt tag;
16baa3c1c6SHong Zhang } Mat_TransMatMultDense;
17baa3c1c6SHong Zhang 
184222ddf1SHong 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;
23faa55883SToby Isaac   PetscInt    alg; /* algorithm used */
24cc48ffa7SToby Isaac } Mat_MatTransMultDense;
25cc48ffa7SToby Isaac 
268965ea79SLois Curfman McInnes typedef struct {
27bad5519cSBarry Smith   Mat         A;                        /* local submatrix */
2813f74950SBarry Smith   PetscMPIInt size;                     /* size of communicator */
2913f74950SBarry Smith   PetscMPIInt rank;                     /* rank of proc in communicator */
30320f2790SHong Zhang 
318965ea79SLois Curfman McInnes   /* The following variables are used for matrix assembly */
32ace3abfcSBarry Smith   PetscBool   donotstash;               /* Flag indicationg if values should be stashed */
338965ea79SLois Curfman McInnes   MPI_Request *send_waits;              /* array of send requests */
348965ea79SLois Curfman McInnes   MPI_Request *recv_waits;              /* array of receive requests */
3513f74950SBarry Smith   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
36ea709b57SSatish Balay   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
3713f74950SBarry Smith   PetscInt    rmax;                     /* maximum message length */
388965ea79SLois Curfman McInnes 
398965ea79SLois Curfman McInnes   /* The following variables are used for matrix-vector products */
408965ea79SLois Curfman McInnes   Vec        lvec;                      /* local vector */
41637a0070SStefano Zampini   PetscSF    Mvctx;                     /* for mat-mult communications */
42ace3abfcSBarry Smith   PetscBool  roworiented;               /* if true, row oriented input (default) */
438949adfdSHong Zhang 
44*5ea7661aSPierre Jolivet   /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
45*5ea7661aSPierre Jolivet   Mat               cmat;      /* matrix representation of a given subset of columns */
466947451fSStefano Zampini   Vec               cvec;      /* vector representation of a given column */
476947451fSStefano Zampini   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
486947451fSStefano Zampini   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
49*5ea7661aSPierre Jolivet   PetscInt          matinuse;  /* if cmat is in use (cbegin = matinuse-1) */
508965ea79SLois Curfman McInnes } Mat_MPIDense;
5197304618SKris Buschelman 
525a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
537dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
544222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat);
55320f2790SHong Zhang 
564222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat);
574222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat);
584222ddf1SHong Zhang 
595242a7b1SHong Zhang #if defined(PETSC_HAVE_ELEMENTAL)
604222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat,Mat,PetscReal,Mat);
614222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat,Mat,Mat);
625242a7b1SHong Zhang #endif
63