Actual source code: mpiaij.h

petsc-master 2020-12-02
Report Typos and Errors

  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;                /* scatter context for vector */
 58:   PetscBool  roworiented;          /* if true, row-oriented input, default true */

 60:   /* The following variables are for MatGetRow() */
 61:   PetscInt    *rowindices;         /* column indices for row */
 62:   PetscScalar *rowvalues;          /* nonzero values in row */
 63:   PetscBool   getrowactive;        /* indicates MatGetRow(), not restored */

 65:   /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation  and nonzero pattern */
 66:   PetscInt *ld;                    /* number of entries per row left of diagona block */

 68:   /* Used by MPICUSPARSE classes */
 69:   void * spptr;

 71: } Mat_MPIAIJ;

 73: PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat);

 75: PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPIAIJ(Mat,MatAssemblyType);

 77: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
 78: PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat);
 79: PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*);
 80: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
 81: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt);
 82: PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring);
 83: PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring);
 84: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 85: PETSC_INTERN PetscErrorCode MatCreateSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 86: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_All(Mat,MatCreateSubMatrixOption,MatReuse,Mat *[]);
 87: PETSC_INTERN PetscErrorCode MatView_MPIAIJ(Mat,PetscViewer);

 89: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*);
 90: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_nonscalable(Mat,IS,IS,PetscInt,MatReuse,Mat*);
 91: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowDist(Mat,IS,IS,IS,MatReuse,Mat*);
 92: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIAIJ_SameRowColDist(Mat,IS,IS,MatReuse,Mat*);
 93: PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*);

 95: PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer);
 96: PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ_Binary(Mat,PetscViewer);
 97: PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat);

 99: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat);
100: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat);

102: PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ(Mat);

104: PETSC_INTERN PetscErrorCode MatProductSymbolic_RARt_MPIAIJ_MPIAIJ(Mat);
105: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt_MPIAIJ_MPIAIJ(Mat);

107: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat);
108: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat,Mat,PetscReal,Mat);
109: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
110: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
111: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);

113: PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat);
114: PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat);

116: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
117: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);

119: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat,Mat,PetscReal,Mat);
120: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat,Mat,PetscReal,Mat);
121: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat,Mat,PetscReal,Mat);
122: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat,Mat,Mat);
123: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat,Mat,Mat);
124: PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat,Mat,Mat);

126: #if defined(PETSC_HAVE_HYPRE)
127: PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
128: #endif
129: #if defined(PETSC_HAVE_SCALAPACK)
130: PETSC_INTERN PetscErrorCode MatConvert_AIJ_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
131: #endif

133: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat);
134: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(void*);
135: PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void*);

137: PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*);
138: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
139: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat(Mat,const PetscInt[],const PetscInt[],const PetscScalar[]);
140: PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Mat,const PetscInt[],const PetscInt[]);
141: PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool);

143: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat);
144: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat);
145: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
146: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat);
147: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat);
148: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
149: PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*);

151: PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat);
152: PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

154: #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16)
155: PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*);
156: #endif
157: PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec);
158: PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*);

160: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*);

162: extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*);
163: extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);

165: PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*);
166: PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat);

168: /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */
169: #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \
170: {\
171:   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ;\
172:   PetscScalar *_aa,_valtmp,*_pa;\
173:   _apJ = apj + api[i];\
174:   /* diagonal portion of A */\
175:   _ai  = ad->i;\
176:   _anz = _ai[i+1] - _ai[i];\
177:   _aj  = ad->j + _ai[i];\
178:   _aa  = ad->a + _ai[i];\
179:   for (_j=0; _j<_anz; _j++) {\
180:     _row = _aj[_j]; \
181:     _pi  = p_loc->i;                             \
182:     _pnz = _pi[_row+1] - _pi[_row];              \
183:     _pj  = p_loc->j + _pi[_row];                 \
184:     _pa  = p_loc->a + _pi[_row];                 \
185:     /* perform sparse axpy */                    \
186:     _valtmp = _aa[_j];                           \
187:     _nextp  = 0; \
188:     for (_k=0; _nextp<_pnz; _k++) {                    \
189:       if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
190:         apa[_k] += _valtmp*_pa[_nextp++];                             \
191:       } \
192:     }                                           \
193:     (void)PetscLogFlops(2.0*_pnz);              \
194:   }                                             \
195:   /* off-diagonal portion of A */               \
196:   if (p_oth){ \
197:     _ai  = ao->i;\
198:     _anz = _ai[i+1] - _ai[i];                   \
199:     _aj  = ao->j + _ai[i];                      \
200:     _aa  = ao->a + _ai[i];                      \
201:     for (_j=0; _j<_anz; _j++) {                 \
202:       _row = _aj[_j];    \
203:       _pi  = p_oth->i;                         \
204:       _pnz = _pi[_row+1] - _pi[_row];          \
205:       _pj  = p_oth->j + _pi[_row];             \
206:       _pa  = p_oth->a + _pi[_row];             \
207:       /* perform sparse axpy */                \
208:       _valtmp = _aa[_j];                       \
209:       _nextp  = 0; \
210:       for (_k=0; _nextp<_pnz; _k++) {          \
211:         if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\
212:           apa[_k] += _valtmp*_pa[_nextp++];    \
213:         }                                      \
214:       }                                        \
215:       (void)PetscLogFlops(2.0*_pnz);           \
216:     } \
217:   }\
218: }

220: #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \
221: {\
222:   PetscInt    _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj;\
223:   PetscScalar *_aa,_valtmp,*_pa;                       \
224:   /* diagonal portion of A */\
225:   _ai  = ad->i;\
226:   _anz = _ai[i+1] - _ai[i];\
227:   _aj  = ad->j + _ai[i];\
228:   _aa  = ad->a + _ai[i];\
229:   for (_j=0; _j<_anz; _j++) {\
230:     _row = _aj[_j]; \
231:     _pi  = p_loc->i;                        \
232:     _pnz = _pi[_row+1] - _pi[_row];         \
233:     _pj  = p_loc->j + _pi[_row];            \
234:     _pa  = p_loc->a + _pi[_row];            \
235:     /* perform dense axpy */                \
236:     _valtmp = _aa[_j];                      \
237:     for (_k=0; _k<_pnz; _k++) {             \
238:       apa[_pj[_k]] += _valtmp*_pa[_k];      \
239:     }                                       \
240:     (void)PetscLogFlops(2.0*_pnz);          \
241:   }                                         \
242:   /* off-diagonal portion of A */           \
243:   if (p_oth){ \
244:     _ai  = ao->i;\
245:     _anz = _ai[i+1] - _ai[i];               \
246:     _aj  = ao->j + _ai[i];                  \
247:     _aa  = ao->a + _ai[i];                  \
248:     for (_j=0; _j<_anz; _j++) {             \
249:       _row = _aj[_j];    \
250:       _pi  = p_oth->i;                      \
251:       _pnz = _pi[_row+1] - _pi[_row];       \
252:       _pj  = p_oth->j + _pi[_row];          \
253:       _pa  = p_oth->a + _pi[_row];          \
254:       /* perform dense axpy */              \
255:       _valtmp = _aa[_j];                    \
256:       for (_k=0; _k<_pnz; _k++) {           \
257:         apa[_pj[_k]] += _valtmp*_pa[_k];    \
258:       }                                     \
259:       (void)PetscLogFlops(2.0*_pnz);        \
260:     }                                       \
261:   }\
262: }

264: #endif