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