Actual source code: mpiaij.h
petsc-3.14.6 2021-03-30
4: #include <../src/mat/impls/aij/seq/aij.h>
6: typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */
7: PetscLayout rowmap;
8: PetscInt **buf_ri,**buf_rj;
9: PetscMPIInt *len_s,*len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */
10: PetscMPIInt nsend,nrecv;
11: PetscInt *bi,*bj; /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */
12: PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */
13: } Mat_Merge_SeqsToMPI;
15: typedef struct { /* used by MatPtAPXXX_MPIAIJ_MPIAIJ() and MatMatMultXXX_MPIAIJ_MPIAIJ() */
16: PetscInt *startsj_s,*startsj_r; /* used by MatGetBrowsOfAoCols_MPIAIJ */
17: PetscScalar *bufa; /* used by MatGetBrowsOfAoCols_MPIAIJ */
18: Mat P_loc,P_oth; /* partial B_seq -- intend to replace B_seq */
19: PetscInt *api,*apj; /* symbolic i and j arrays of the local product A_loc*B_seq */
20: PetscScalar *apv;
21: MatReuse reuse; /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */
22: PetscScalar *apa; /* tmp array for store a row of A*P used in MatMatMult() */
23: Mat A_loc; /* used by MatTransposeMatMult(), contains api and apj */
24: ISLocalToGlobalMapping ltog; /* mapping from local column indices to global column indices for A_loc */
25: Mat Pt; /* used by MatTransposeMatMult(), Pt = P^T */
26: Mat Rd,Ro,AP_loc,C_loc,C_oth;
27: PetscInt algType; /* implementation algorithm */
28: PetscSF sf; /* use it to communicate remote part of C */
29: PetscInt *c_othi,*c_rmti;
31: Mat_Merge_SeqsToMPI *merge;
32: } Mat_APMPI;
34: typedef struct {
35: Mat A,B; /* local submatrices: A (diag part),
36: B (off-diag part) */
37: PetscMPIInt size; /* size of communicator */
38: PetscMPIInt rank; /* rank of proc in communicator */
40: /* The following variables are used for matrix assembly */
41: PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */
42: MPI_Request *send_waits; /* array of send requests */
43: MPI_Request *recv_waits; /* array of receive requests */
44: PetscInt nsends,nrecvs; /* numbers of sends and receives */
45: PetscScalar *svalues,*rvalues; /* sending and receiving data */
46: PetscInt rmax; /* maximum message length */
47: #if defined(PETSC_USE_CTABLE)
48: PetscTable colmap;
49: #else
50: PetscInt *colmap; /* local col number of off-diag col */
51: #endif
52: PetscInt *garray; /* global index of all off-processor columns */
54: /* The following variables are used for matrix-vector products */
55: Vec lvec; /* local vector */
56: Vec diag;
57: VecScatter Mvctx,Mvctx_mpi1; /* scatter context for vector */
58: PetscBool Mvctx_mpi1_flg; /* if true, additional Mvctx_mpi1 is requested for mat-mat ops, default false */
59: PetscBool roworiented; /* if true, row-oriented input, default true */
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 */
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 */
69: /* Used by MPICUSPARSE classes */
70: void * spptr;
72: } Mat_MPIAIJ;
74: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);
76: PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType);
78: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
79: PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
80: PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
81: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
82: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt);
83: PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
84: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
85: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
86: PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
87: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat,MatCreateSubMatrixOption,MatReuse,Mat *[]);
88: PETSC_INTERN PetscErrorCode MatView_MPIAIJ(Mat,PetscViewer);
90: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
91: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat,IS,IS,PetscInt,MatReuse,Mat*);
92: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat,IS,IS,IS,MatReuse,Mat*);
93: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat,IS,IS,MatReuse,Mat*);
94: PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);
96: PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
97: PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ_Binary(Mat,PetscViewer);
98: PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);
100: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat);
101: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat);
103: PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ(Mat);
105: PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_MPIAIJ_MPIAIJ(Mat);
106: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt_MPIAIJ_MPIAIJ(Mat);
108: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat);
109: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat,Mat,PetscReal,Mat);
110: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
111: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
112: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
114: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat);
115: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);
117: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
118: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
120: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat,Mat,PetscReal,Mat);
121: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat,Mat,PetscReal,Mat);
122: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat,Mat,PetscReal,Mat);
123: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat,Mat,Mat);
124: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat,Mat,Mat);
125: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat,Mat,Mat);
127: #if defined(PETSC_HAVE_HYPRE)
128: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
129: #endif
130: #if defined(PETSC_HAVE_SCALAPACK)
131: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
132: #endif
134: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
135: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(void*);
136: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void*);
138: PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
139: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
140: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
141: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat,const PetscInt[],const PetscInt[]);
142: PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);
144: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat);
145: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
146: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
147: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
148: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
149: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
150: PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);
152: PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
153: PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);
155: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
156: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
157: #endif
158: PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
159: PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);
161: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);
163: extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
164: extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
166: PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
167: PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);
169: /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
170: #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
171: {\
172: PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ;\
173: PetscScalar *_aa,_valtmp,*_pa;\
174: _apJ = apj + api[i];\
175: /* diagonal portion of A */\
176: _ai = ad->i;\
177: _anz = _ai[i+1] - _ai[i];\
178: _aj = ad->j + _ai[i];\
179: _aa = ad->a + _ai[i];\
180: for (_j=0; _j<_anz; _j++) {\
181: _row = _aj[_j]; \
182: _pi = p_loc->i; \
183: _pnz = _pi[_row+1] - _pi[_row]; \
184: _pj = p_loc->j + _pi[_row]; \
185: _pa = p_loc->a + _pi[_row]; \
186: /* perform sparse axpy */ \
187: _valtmp = _aa[_j]; \
188: _nextp = 0; \
189: for (_k=0; _nextp<_pnz; _k++) { \
190: if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
191: apa[_k] += _valtmp*_pa[_nextp++]; \
192: } \
193: } \
194: (void)PetscLogFlops(2.0*_pnz); \
195: } \
196: /* off-diagonal portion of A */ \
197: if (p_oth){ \
198: _ai = ao->i;\
199: _anz = _ai[i+1] - _ai[i]; \
200: _aj = ao->j + _ai[i]; \
201: _aa = ao->a + _ai[i]; \
202: for (_j=0; _j<_anz; _j++) { \
203: _row = _aj[_j]; \
204: _pi = p_oth->i; \
205: _pnz = _pi[_row+1] - _pi[_row]; \
206: _pj = p_oth->j + _pi[_row]; \
207: _pa = p_oth->a + _pi[_row]; \
208: /* perform sparse axpy */ \
209: _valtmp = _aa[_j]; \
210: _nextp = 0; \
211: for (_k=0; _nextp<_pnz; _k++) { \
212: if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
213: apa[_k] += _valtmp*_pa[_nextp++]; \
214: } \
215: } \
216: (void)PetscLogFlops(2.0*_pnz); \
217: } \
218: }\
219: }
221: #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
222: {\
223: PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;\
224: PetscScalar *_aa,_valtmp,*_pa; \
225: /* diagonal portion of A */\
226: _ai = ad->i;\
227: _anz = _ai[i+1] - _ai[i];\
228: _aj = ad->j + _ai[i];\
229: _aa = ad->a + _ai[i];\
230: for (_j=0; _j<_anz; _j++) {\
231: _row = _aj[_j]; \
232: _pi = p_loc->i; \
233: _pnz = _pi[_row+1] - _pi[_row]; \
234: _pj = p_loc->j + _pi[_row]; \
235: _pa = p_loc->a + _pi[_row]; \
236: /* perform dense axpy */ \
237: _valtmp = _aa[_j]; \
238: for (_k=0; _k<_pnz; _k++) { \
239: apa[_pj[_k]] += _valtmp*_pa[_k]; \
240: } \
241: (void)PetscLogFlops(2.0*_pnz); \
242: } \
243: /* off-diagonal portion of A */ \
244: if (p_oth){ \
245: _ai = ao->i;\
246: _anz = _ai[i+1] - _ai[i]; \
247: _aj = ao->j + _ai[i]; \
248: _aa = ao->a + _ai[i]; \
249: for (_j=0; _j<_anz; _j++) { \
250: _row = _aj[_j]; \
251: _pi = p_oth->i; \
252: _pnz = _pi[_row+1] - _pi[_row]; \
253: _pj = p_oth->j + _pi[_row]; \
254: _pa = p_oth->a + _pi[_row]; \
255: /* perform dense axpy */ \
256: _valtmp = _aa[_j]; \
257: for (_k=0; _k<_pnz; _k++) { \
258: apa[_pj[_k]] += _valtmp*_pa[_k]; \
259: } \
260: (void)PetscLogFlops(2.0*_pnz); \
261: } \
262: }\
263: }
265: #endif