xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
1 
2 #if !defined(__MPIAIJ_H)
3 #define __MPIAIJ_H
4 
5 #include <../src/mat/impls/aij/seq/aij.h>
6 #include <petscctable.h>
7 
8 typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
9   PetscLayout rowmap;
10   PetscInt    **buf_ri,**buf_rj;
11   PetscMPIInt *len_s,*len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
12   PetscMPIInt nsend,nrecv;
13   PetscInt    *bi,*bj;    /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
14   PetscInt    *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
15   PetscErrorCode (*destroy)(Mat);
16   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
17 } Mat_Merge_SeqsToMPI;
18 
19 typedef struct { /* used by MatPtAP_MPIAIJ_MPIAIJ() and MatMatMult_MPIAIJ_MPIAIJ() */
20   PetscInt    *startsj_s,*startsj_r;    /* used by MatGetBrowsOfAoCols_MPIAIJ */
21   PetscScalar *bufa;                    /* used by MatGetBrowsOfAoCols_MPIAIJ */
22   Mat         P_loc,P_oth;     /* partial B_seq -- intend to replace B_seq */
23   PetscInt    *api,*apj;       /* symbolic i and j arrays of the local product A_loc*B_seq */
24   PetscInt    rmax;            /* max num of nnz in a local row of the matrix product */
25   MatReuse    reuse;           /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
26   PetscScalar *apa;            /* tmp array for store a row of A*P used in MatMatMult() */
27   Mat         A_loc;           /* used by MatTransposeMatMult(), contains api and apj */
28   Mat         Pt;              /* used by MatTransposeMatMult(), Pt = P^T */
29   PetscBool   scalable;        /* flag determines scalable or non-scalable implementation */
30 
31   Mat_Merge_SeqsToMPI *merge;
32   PetscErrorCode (*destroy)(Mat);
33   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
34 } Mat_PtAPMPI;
35 
36 typedef struct {
37   Mat A,B;                             /* local submatrices: A (diag part),
38                                            B (off-diag part) */
39   PetscMPIInt size;                     /* size of communicator */
40   PetscMPIInt rank;                     /* rank of proc in communicator */
41 
42   /* The following variables are used for matrix assembly */
43   PetscBool   donotstash;               /* PETSC_TRUE if off processor entries dropped */
44   MPI_Request *send_waits;              /* array of send requests */
45   MPI_Request *recv_waits;              /* array of receive requests */
46   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
47   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
48   PetscInt    rmax;                     /* maximum message length */
49 #if defined(PETSC_USE_CTABLE)
50   PetscTable colmap;
51 #else
52   PetscInt *colmap;                     /* local col number of off-diag col */
53 #endif
54   PetscInt *garray;                     /* global index of all off-processor columns */
55 
56   /* The following variables are used for matrix-vector products */
57   Vec        lvec;                 /* local vector */
58   Vec        diag;
59   VecScatter Mvctx;                /* scatter context for vector */
60   PetscBool  roworiented;          /* if true, row-oriented input, default true */
61 
62   /* The following variables are for MatGetRow() */
63   PetscInt    *rowindices;         /* column indices for row */
64   PetscScalar *rowvalues;          /* nonzero values in row */
65   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */
66 
67   /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation  and nonzero pattern */
68   PetscInt *ld;                    /* number of entries per row left of diagona block */
69 
70   /* Used by MatMatMult() and MatPtAP() */
71   Mat_PtAPMPI *ptap;
72 
73   /* used by MatMatMatMult() */
74   Mat_MatMatMatMult *matmatmatmult;
75 
76   /* Used by MPICUSP and MPICUSPARSE classes */
77   void * spptr;
78 
79 } Mat_MPIAIJ;
80 
81 PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
82 
83 PETSC_INTERN PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring);
84 PETSC_INTERN PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*);
85 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
86 PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
87 PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
88 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
89 PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
90 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
91 PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
92 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat,MatGetSubMatrixOption,MatReuse,Mat *[]);
93 PETSC_INTERN PetscErrorCode MatGetSubMatricesParallel_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
94 
95 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
96 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_Private (Mat,IS,IS,PetscInt,MatReuse,Mat*);
97 PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
98 
99 PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
100 PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
101 PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
102 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
103 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
104 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
105 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
106 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
107 
108 PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
109 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*);
110 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
111 
112 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
113 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
114 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
115 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
116 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
117 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat);
118 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
119 
120 PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
121 PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
122 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat);
123 PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*);
124 PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
125 
126 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
127 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
128 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
129 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
130 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
131 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
132 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
133 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
134 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
135 PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
136 
137 extern PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
138 
139 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128)
140 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
141 #endif
142 PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
143 PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
144 
145 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
146 
147 extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
148 extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
149 
150 #endif
151