xref: /petsc/src/mat/impls/dense/mpi/mpidense.h (revision 8a86b0881472c00bf0706e4718310bd7d50b8cec)
1 
2 #include <../src/mat/impls/dense/seq/dense.h>
3 #include <petscsf.h>
4 
5 /*  Data stuctures for basic parallel dense matrix  */
6 
7 typedef struct { /* used by MatMatMultxxx_MPIDense_MPIDense() */
8   Mat            Ae,Be,Ce;           /* matrix in Elemental format */
9   PetscErrorCode (*destroy)(Mat);
10 } Mat_MatMultDense;
11 
12 typedef struct { /* used by MatTransposeMatMultXXX_MPIDense_MPIDense() */
13   PetscScalar    *sendbuf;
14   Mat            atb;
15   PetscMPIInt    *recvcounts;
16   PetscErrorCode (*destroy)(Mat);
17   PetscMPIInt    tag;
18 } Mat_TransMatMultDense;
19 
20 typedef struct { /* used by MatMatTransposeMultxxx_MPIDense_MPIDense() */
21   PetscScalar    *buf[2];
22   PetscMPIInt    tag;
23   PetscMPIInt    *recvcounts;
24   PetscMPIInt    *recvdispls;
25   PetscErrorCode (*destroy)(Mat);
26   PetscInt       alg; /* algorithm used */
27 } Mat_MatTransMultDense;
28 
29 typedef struct {
30   Mat         A;                        /* local submatrix */
31   PetscMPIInt size;                     /* size of communicator */
32   PetscMPIInt rank;                     /* rank of proc in communicator */
33 
34   /* The following variables are used for matrix assembly */
35   PetscBool   donotstash;               /* Flag indicationg if values should be stashed */
36   MPI_Request *send_waits;              /* array of send requests */
37   MPI_Request *recv_waits;              /* array of receive requests */
38   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
39   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
40   PetscInt    rmax;                     /* maximum message length */
41 
42   /* The following variables are used for matrix-vector products */
43   Vec        lvec;                      /* local vector */
44   PetscSF    Mvctx;                     /* for mat-mult communications */
45   PetscBool  roworiented;               /* if true, row oriented input (default) */
46 
47   /* Support for MatDenseGetColumnVec */
48   Vec               cvec;      /* vector representation of a given column */
49   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
50   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
51 
52   Mat_MatTransMatMult   *atb;           /* used by MatTransposeMatMultxxx_MPIAIJ_MPIDense */
53   Mat_TransMatMultDense *atbdense;      /* used by MatTransposeMatMultxxx_MPIDense_MPIDense */
54   Mat_MatMultDense      *abdense;       /* used by MatMatMultxxx_MPIDense_MPIDense */
55   Mat_MatTransMultDense *abtdense;      /* used by MatMatTransposeMultxxx_MPIDense_MPIDense */
56 } Mat_MPIDense;
57 
58 PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer);
59 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
60 PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
61 PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*);
62 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
63 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat);
64 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
65 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
66 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat);
67 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
68 
69 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat);
70 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat);
71 
72 #if defined(PETSC_HAVE_ELEMENTAL)
73 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
74 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat,Mat,PetscReal,Mat);
75 PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat,Mat,Mat);
76 #endif
77