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