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 PetscInt rmax; /* max num of nnz in a local row of the matrix product */ 25 MatReuse reuse; /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 26 PetscScalar *apa; /* tmp array for store a row of A*P used in MatMatMult() */ 27 Mat A_loc; /* used by MatTransposeMatMult(), contains api and apj */ 28 Mat Pt; /* used by MatTransposeMatMult(), Pt = P^T */ 29 PetscBool scalable; /* flag determines scalable or non-scalable implementation */ 30 Mat Rd,Ro,AP_loc,C_loc,C_oth; 31 32 Mat_Merge_SeqsToMPI *merge; 33 PetscErrorCode (*destroy)(Mat); 34 PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 35 } Mat_PtAPMPI; 36 37 typedef struct { 38 Mat A,B; /* local submatrices: A (diag part), 39 B (off-diag part) */ 40 PetscMPIInt size; /* size of communicator */ 41 PetscMPIInt rank; /* rank of proc in communicator */ 42 43 /* The following variables are used for matrix assembly */ 44 PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */ 45 MPI_Request *send_waits; /* array of send requests */ 46 MPI_Request *recv_waits; /* array of receive requests */ 47 PetscInt nsends,nrecvs; /* numbers of sends and receives */ 48 PetscScalar *svalues,*rvalues; /* sending and receiving data */ 49 PetscInt rmax; /* maximum message length */ 50 #if defined(PETSC_USE_CTABLE) 51 PetscTable colmap; 52 #else 53 PetscInt *colmap; /* local col number of off-diag col */ 54 #endif 55 PetscInt *garray; /* global index of all off-processor columns */ 56 57 /* The following variables are used for matrix-vector products */ 58 Vec lvec; /* local vector */ 59 Vec diag; 60 VecScatter Mvctx; /* scatter context for vector */ 61 PetscBool roworiented; /* if true, row-oriented input, default true */ 62 63 /* The following variables are for MatGetRow() */ 64 PetscInt *rowindices; /* column indices for row */ 65 PetscScalar *rowvalues; /* nonzero values in row */ 66 PetscBool getrowactive; /* indicates MatGetRow(), not restored */ 67 68 /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation and nonzero pattern */ 69 PetscInt *ld; /* number of entries per row left of diagona block */ 70 71 /* Used by MatMatMult() and MatPtAP() */ 72 Mat_PtAPMPI *ptap; 73 74 /* used by MatMatMatMult() */ 75 Mat_MatMatMatMult *matmatmatmult; 76 77 /* Used by MPICUSP and MPICUSPARSE classes */ 78 void * spptr; 79 80 } Mat_MPIAIJ; 81 82 PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat); 83 84 PETSC_INTERN PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring); 85 PETSC_INTERN PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*); 86 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat); 87 PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat); 88 PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*); 89 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt); 90 PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ_Scalable(Mat,PetscInt,IS [],PetscInt); 91 PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring); 92 PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring); 93 PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 94 PETSC_INTERN PetscErrorCode MatGetSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 95 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat,MatGetSubMatrixOption,MatReuse,Mat *[]); 96 97 98 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*); 99 PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_Private (Mat,IS,IS,PetscInt,MatReuse,Mat*); 100 PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*); 101 102 PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer); 103 PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); 104 PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 105 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 106 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 107 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 108 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 109 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 110 111 PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 112 PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*); 113 PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat); 114 115 PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 116 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 117 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 118 119 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat,Mat,PetscReal,Mat*); 120 PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat,Mat,Mat); 121 122 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat); 123 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 124 125 PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 126 PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 127 PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 128 PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*); 129 PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool); 130 131 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 132 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 133 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 134 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 135 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 136 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat); 137 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 138 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 139 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 140 PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*); 141 142 PETSC_INTERN PetscErrorCode MatSetFromOptions_MPIAIJ(PetscOptions*,Mat); 143 PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 144 145 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) 146 PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*); 147 #endif 148 PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec); 149 PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*); 150 151 PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 152 153 extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*); 154 extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 155 156 PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*); 157 PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat); 158 159 /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth */ 160 #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \ 161 {\ 162 PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj; \ 163 PetscScalar *_aa,_valtmp,*_pa; \ 164 /* diagonal portion of A */\ 165 _ai = ad->i;\ 166 _anz = _ai[i+1] - _ai[i];\ 167 _aj = ad->j + _ai[i];\ 168 _aa = ad->a + _ai[i];\ 169 for (_j=0; _j<_anz; _j++) {\ 170 _row = _aj[_j]; \ 171 _pi = p_loc->i; \ 172 _pnz = _pi[_row+1] - _pi[_row]; \ 173 _pj = p_loc->j + _pi[_row]; \ 174 _pa = p_loc->a + _pi[_row]; \ 175 /* perform dense axpy */ \ 176 _valtmp = _aa[_j]; \ 177 for (_k=0; _k<_pnz; _k++) { \ 178 apa[_pj[_k]] += _valtmp*_pa[_k]; \ 179 } \ 180 PetscLogFlops(2.0*_pnz); \ 181 } \ 182 /* off-diagonal portion of A */ \ 183 _ai = ao->i;\ 184 _anz = _ai[i+1] - _ai[i]; \ 185 _aj = ao->j + _ai[i]; \ 186 _aa = ao->a + _ai[i]; \ 187 for (_j=0; _j<_anz; _j++) { \ 188 _row = _aj[_j]; \ 189 _pi = p_oth->i; \ 190 _pnz = _pi[_row+1] - _pi[_row]; \ 191 _pj = p_oth->j + _pi[_row]; \ 192 _pa = p_oth->a + _pi[_row]; \ 193 /* perform dense axpy */ \ 194 _valtmp = _aa[_j]; \ 195 for (_k=0; _k<_pnz; _k++) { \ 196 apa[_pj[_k]] += _valtmp*_pa[_k]; \ 197 } \ 198 PetscLogFlops(2.0*_pnz); \ 199 } \ 200 } 201 #endif 202