Actual source code: mpiaij.h
petsc-3.11.4 2019-09-28
5: #include <../src/mat/impls/aij/seq/aij.h>
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;
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: PetscScalar *apv;
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: ISLocalToGlobalMapping ltog; /* mapping from local column indices to global column indices for A_loc */
28: Mat Pt; /* used by MatTransposeMatMult(), Pt = P^T */
29: PetscBool freestruct; /* flag for MatFreeIntermediateDataStructures() */
30: Mat Rd,Ro,AP_loc,C_loc,C_oth;
31: PetscInt algType; /* implementation algorithm */
33: Mat_Merge_SeqsToMPI *merge;
34: PetscErrorCode (*destroy)(Mat);
35: PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
36: PetscErrorCode (*view)(Mat,PetscViewer);
37: } Mat_APMPI;
39: typedef struct {
40: Mat A,B; /* local submatrices: A (diag part),
41: B (off-diag part) */
42: PetscMPIInt size; /* size of communicator */
43: PetscMPIInt rank; /* rank of proc in communicator */
45: /* The following variables are used for matrix assembly */
46: PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */
47: MPI_Request *send_waits; /* array of send requests */
48: MPI_Request *recv_waits; /* array of receive requests */
49: PetscInt nsends,nrecvs; /* numbers of sends and receives */
50: PetscScalar *svalues,*rvalues; /* sending and receiving data */
51: PetscInt rmax; /* maximum message length */
52: #if defined(PETSC_USE_CTABLE)
53: PetscTable colmap;
54: #else
55: PetscInt *colmap; /* local col number of off-diag col */
56: #endif
57: PetscInt *garray; /* global index of all off-processor columns */
59: /* The following variables are used for matrix-vector products */
60: Vec lvec; /* local vector */
61: Vec diag;
62: VecScatter Mvctx,Mvctx_mpi1; /* scatter context for vector */
63: PetscBool Mvctx_mpi1_flg; /* if true, additional Mvctx_mpi1 is requested for mat-mat ops, default false */
64: PetscBool roworiented; /* if true, row-oriented input, default true */
66: /* The following variables are for MatGetRow() */
67: PetscInt *rowindices; /* column indices for row */
68: PetscScalar *rowvalues; /* nonzero values in row */
69: PetscBool getrowactive; /* indicates MatGetRow(), not restored */
71: /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation and nonzero pattern */
72: PetscInt *ld; /* number of entries per row left of diagona block */
74: /* Used by MatMatMult() and MatPtAP() */
75: Mat_APMPI *ap;
77: /* used by MatMatMatMult() */
78: Mat_MatMatMatMult *matmatmatmult;
80: /* Used by MPICUSP and MPICUSPARSE classes */
81: void * spptr;
83: } Mat_MPIAIJ;
85: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
87: PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType);
89: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
90: PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
91: PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
92: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
93: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt);
94: PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
95: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
96: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
97: PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
98: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat,MatCreateSubMatrixOption,MatReuse,Mat *[]);
99: PETSC_INTERN PetscErrorCode MatView_MPIAIJ(Mat,PetscViewer);
101: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
102: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat,IS,IS,PetscInt,MatReuse,Mat*);
103: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat,IS,IS,IS,MatReuse,Mat*);
104: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat,IS,IS,MatReuse,Mat*);
105: PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
107: PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
108: PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ_Binary(Mat,PetscViewer);
109: PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
110: PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
111: PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
112: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
113: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat,Mat,PetscReal,Mat*);
114: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
115: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
116: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
118: PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
119: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*);
120: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
122: PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
123: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
124: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
126: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat,Mat,PetscReal,Mat*);
127: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat,Mat,Mat);
128: PETSC_INTERN PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_AP(Mat);
129: PETSC_INTERN PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_BC(Mat);
131: #if defined(PETSC_HAVE_HYPRE)
132: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
133: #endif
135: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat);
136: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
138: PETSC_INTERN PetscErrorCode MatRARt_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
140: PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
141: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
142: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
143: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat,const PetscInt[],const PetscInt[]);
144: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat);
145: PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*);
146: PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
148: PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
149: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*);
150: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
151: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
152: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
153: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
154: PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
155: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
156: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
157: PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
159: PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
160: PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
162: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
163: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
164: #endif
165: PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
166: PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
168: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
170: extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
171: extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
173: PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
174: PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);
176: /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
177: #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
178: {\
179: PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ; \
180: PetscScalar *_aa,_valtmp,*_pa; \
181: _apJ = apj + api[i];\
182: /* diagonal portion of A */\
183: _ai = ad->i;\
184: _anz = _ai[i+1] - _ai[i];\
185: _aj = ad->j + _ai[i];\
186: _aa = ad->a + _ai[i];\
187: for (_j=0; _j<_anz; _j++) {\
188: _row = _aj[_j]; \
189: _pi = p_loc->i; \
190: _pnz = _pi[_row+1] - _pi[_row]; \
191: _pj = p_loc->j + _pi[_row]; \
192: _pa = p_loc->a + _pi[_row]; \
193: /* perform sparse axpy */ \
194: _valtmp = _aa[_j]; \
195: _nextp = 0; \
196: for (_k=0; _nextp<_pnz; _k++) { \
197: if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \
198: apa[_k] += _valtmp*_pa[_nextp++]; \
199: } \
200: } \
201: (void)PetscLogFlops(2.0*_pnz); \
202: } \
203: /* off-diagonal portion of A */ \
204: _ai = ao->i;\
205: _anz = _ai[i+1] - _ai[i]; \
206: _aj = ao->j + _ai[i]; \
207: _aa = ao->a + _ai[i]; \
208: for (_j=0; _j<_anz; _j++) { \
209: _row = _aj[_j]; \
210: _pi = p_oth->i; \
211: _pnz = _pi[_row+1] - _pi[_row]; \
212: _pj = p_oth->j + _pi[_row]; \
213: _pa = p_oth->a + _pi[_row]; \
214: /* perform sparse axpy */ \
215: _valtmp = _aa[_j]; \
216: _nextp = 0; \
217: for (_k=0; _nextp<_pnz; _k++) { \
218: if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
219: apa[_k] += _valtmp*_pa[_nextp++]; \
220: } \
221: } \
222: (void)PetscLogFlops(2.0*_pnz); \
223: } \
224: }
226: #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
227: {\
228: PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj; \
229: PetscScalar *_aa,_valtmp,*_pa; \
230: /* diagonal portion of A */\
231: _ai = ad->i;\
232: _anz = _ai[i+1] - _ai[i];\
233: _aj = ad->j + _ai[i];\
234: _aa = ad->a + _ai[i];\
235: for (_j=0; _j<_anz; _j++) {\
236: _row = _aj[_j]; \
237: _pi = p_loc->i; \
238: _pnz = _pi[_row+1] - _pi[_row]; \
239: _pj = p_loc->j + _pi[_row]; \
240: _pa = p_loc->a + _pi[_row]; \
241: /* perform dense axpy */ \
242: _valtmp = _aa[_j]; \
243: for (_k=0; _k<_pnz; _k++) { \
244: apa[_pj[_k]] += _valtmp*_pa[_k]; \
245: } \
246: (void)PetscLogFlops(2.0*_pnz); \
247: } \
248: /* off-diagonal portion of A */ \
249: _ai = ao->i;\
250: _anz = _ai[i+1] - _ai[i]; \
251: _aj = ao->j + _ai[i]; \
252: _aa = ao->a + _ai[i]; \
253: for (_j=0; _j<_anz; _j++) { \
254: _row = _aj[_j]; \
255: _pi = p_oth->i; \
256: _pnz = _pi[_row+1] - _pi[_row]; \
257: _pj = p_oth->j + _pi[_row]; \
258: _pa = p_oth->a + _pi[_row]; \
259: /* perform dense axpy */ \
260: _valtmp = _aa[_j]; \
261: for (_k=0; _k<_pnz; _k++) { \
262: apa[_pj[_k]] += _valtmp*_pa[_k]; \
263: } \
264: (void)PetscLogFlops(2.0*_pnz); \
265: } \
266: }
268: #endif