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