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_Private (Mat,IS,IS,PetscInt,MatReuse,Mat*); 97 PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*); 98 99 PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer); 100 PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); 101 PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 102 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 103 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 104 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 105 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 106 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 107 108 PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 109 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*); 110 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat); 111 112 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 113 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 114 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 115 116 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat,Mat,PetscReal,Mat*); 117 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat,Mat,Mat); 118 119 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat); 120 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 121 122 PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 123 PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 124 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 125 PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*); 126 PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool); 127 128 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 129 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 130 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 131 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 132 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 133 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat); 134 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 135 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 136 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 137 PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*); 138 139 PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptionItems*,Mat); 140 PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 141 142 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) && !defined(PETSC_USE_REAL___FP16) 143 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*); 144 #endif 145 PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec); 146 PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*); 147 148 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 149 150 extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*); 151 extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 152 153 PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*); 154 PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat); 155 156 /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth using sparse axpy */ 157 #define AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa) \ 158 {\ 159 PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj,_nextp,*_apJ; \ 160 PetscScalar *_aa,_valtmp,*_pa; \ 161 _apJ = apj + api[i];\ 162 /* diagonal portion of A */\ 163 _ai = ad->i;\ 164 _anz = _ai[i+1] - _ai[i];\ 165 _aj = ad->j + _ai[i];\ 166 _aa = ad->a + _ai[i];\ 167 for (_j=0; _j<_anz; _j++) {\ 168 _row = _aj[_j]; \ 169 _pi = p_loc->i; \ 170 _pnz = _pi[_row+1] - _pi[_row]; \ 171 _pj = p_loc->j + _pi[_row]; \ 172 _pa = p_loc->a + _pi[_row]; \ 173 /* perform sparse axpy */ \ 174 _valtmp = _aa[_j]; \ 175 _nextp = 0; \ 176 for (_k=0; _nextp<_pnz; _k++) { \ 177 if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */ \ 178 apa[_k] += _valtmp*_pa[_nextp++]; \ 179 } \ 180 } \ 181 PetscLogFlops(2.0*_pnz); \ 182 } \ 183 /* off-diagonal portion of A */ \ 184 _ai = ao->i;\ 185 _anz = _ai[i+1] - _ai[i]; \ 186 _aj = ao->j + _ai[i]; \ 187 _aa = ao->a + _ai[i]; \ 188 for (_j=0; _j<_anz; _j++) { \ 189 _row = _aj[_j]; \ 190 _pi = p_oth->i; \ 191 _pnz = _pi[_row+1] - _pi[_row]; \ 192 _pj = p_oth->j + _pi[_row]; \ 193 _pa = p_oth->a + _pi[_row]; \ 194 /* perform sparse axpy */ \ 195 _valtmp = _aa[_j]; \ 196 _nextp = 0; \ 197 for (_k=0; _nextp<_pnz; _k++) { \ 198 if (_apJ[_k] == _pj[_nextp]) { /* column of AP == column of P */\ 199 apa[_k] += _valtmp*_pa[_nextp++]; \ 200 } \ 201 } \ 202 PetscLogFlops(2.0*_pnz); \ 203 } \ 204 } 205 206 #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \ 207 {\ 208 PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj; \ 209 PetscScalar *_aa,_valtmp,*_pa; \ 210 /* diagonal portion of A */\ 211 _ai = ad->i;\ 212 _anz = _ai[i+1] - _ai[i];\ 213 _aj = ad->j + _ai[i];\ 214 _aa = ad->a + _ai[i];\ 215 for (_j=0; _j<_anz; _j++) {\ 216 _row = _aj[_j]; \ 217 _pi = p_loc->i; \ 218 _pnz = _pi[_row+1] - _pi[_row]; \ 219 _pj = p_loc->j + _pi[_row]; \ 220 _pa = p_loc->a + _pi[_row]; \ 221 /* perform dense axpy */ \ 222 _valtmp = _aa[_j]; \ 223 for (_k=0; _k<_pnz; _k++) { \ 224 apa[_pj[_k]] += _valtmp*_pa[_k]; \ 225 } \ 226 PetscLogFlops(2.0*_pnz); \ 227 } \ 228 /* off-diagonal portion of A */ \ 229 _ai = ao->i;\ 230 _anz = _ai[i+1] - _ai[i]; \ 231 _aj = ao->j + _ai[i]; \ 232 _aa = ao->a + _ai[i]; \ 233 for (_j=0; _j<_anz; _j++) { \ 234 _row = _aj[_j]; \ 235 _pi = p_oth->i; \ 236 _pnz = _pi[_row+1] - _pi[_row]; \ 237 _pj = p_oth->j + _pi[_row]; \ 238 _pa = p_oth->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 PetscLogFlops(2.0*_pnz); \ 245 } \ 246 } 247 248 #endif 249