18a729477SBarry Smith 23369ce9aSBarry Smith #if !defined(__MPIAIJ_H) 33369ce9aSBarry Smith #define __MPIAIJ_H 43369ce9aSBarry Smith 5c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 682c86c8fSBarry Smith #include <petscctable.h> 78a729477SBarry Smith 890431a8fSHong Zhang typedef struct { /* used by MatCreateMPIAIJSumSeqAIJ for reusing the merged matrix */ 926283091SBarry Smith PetscLayout rowmap; 10de0260b3SHong Zhang PetscInt **buf_ri,**buf_rj; 110febcb4bSHong Zhang PetscMPIInt *len_s,*len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 120febcb4bSHong Zhang PetscMPIInt nsend,nrecv; 13fc08c53fSHong Zhang PetscInt *bi,*bj; /* i and j array of the local portion of mpi C (matrix product) - rename to ci, cj! */ 14de0260b3SHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 15dce485f0SBarry Smith PetscErrorCode (*destroy)(Mat); 16dce485f0SBarry Smith PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 17b90dcfe3SHong Zhang } Mat_Merge_SeqsToMPI; 18b90dcfe3SHong Zhang 19a1a4e84aSHong Zhang typedef struct { /* used by MatPtAP_MPIAIJ_MPIAIJ() and MatMatMult_MPIAIJ_MPIAIJ() */ 20b7f45c76SHong Zhang PetscInt *startsj_s,*startsj_r; /* used by MatGetBrowsOfAoCols_MPIAIJ */ 21b7f45c76SHong Zhang PetscScalar *bufa; /* used by MatGetBrowsOfAoCols_MPIAIJ */ 22a1a4e84aSHong Zhang Mat P_loc,P_oth; /* partial B_seq -- intend to replace B_seq */ 23a1a4e84aSHong Zhang PetscInt *api,*apj; /* symbolic i and j arrays of the local product A_loc*B_seq */ 24d6ab1dc0SHong Zhang PetscInt rmax; /* max num of nnz in a local row of the matrix product */ 25b7f45c76SHong Zhang MatReuse reuse; /* flag to skip MatGetBrowsOfAoCols_MPIAIJ() and MatMPIAIJGetLocalMat() in 1st call of MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 26598bc09dSHong Zhang PetscScalar *apa; /* tmp array for store a row of A*P used in MatMatMult() */ 27c5af039cSHong Zhang Mat A_loc; /* used by MatTransposeMatMult(), contains api and apj */ 28ab1b013aSHong Zhang Mat Pt; /* used by MatTransposeMatMult(), Pt = P^T */ 29414894bdSHong Zhang PetscBool scalable; /* flag determines scalable or non-scalable implementation */ 300d3441aeSHong Zhang Mat Rd,Ro,AP,AP_loc; 31f8487c73SHong Zhang 32f8487c73SHong Zhang Mat_Merge_SeqsToMPI *merge; 33f8487c73SHong Zhang PetscErrorCode (*destroy)(Mat); 344ae0847bSHong Zhang PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 35f8487c73SHong Zhang } Mat_PtAPMPI; 36f8487c73SHong Zhang 37f8487c73SHong Zhang typedef struct { 38f8487c73SHong Zhang Mat A,B; /* local submatrices: A (diag part), 39f8487c73SHong Zhang B (off-diag part) */ 40f8487c73SHong Zhang PetscMPIInt size; /* size of communicator */ 41f8487c73SHong Zhang PetscMPIInt rank; /* rank of proc in communicator */ 42f8487c73SHong Zhang 43f8487c73SHong Zhang /* The following variables are used for matrix assembly */ 44f8487c73SHong Zhang PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */ 45f8487c73SHong Zhang MPI_Request *send_waits; /* array of send requests */ 46f8487c73SHong Zhang MPI_Request *recv_waits; /* array of receive requests */ 47f8487c73SHong Zhang PetscInt nsends,nrecvs; /* numbers of sends and receives */ 48f8487c73SHong Zhang PetscScalar *svalues,*rvalues; /* sending and receiving data */ 49f8487c73SHong Zhang PetscInt rmax; /* maximum message length */ 50f8487c73SHong Zhang #if defined(PETSC_USE_CTABLE) 51f8487c73SHong Zhang PetscTable colmap; 52f8487c73SHong Zhang #else 53f8487c73SHong Zhang PetscInt *colmap; /* local col number of off-diag col */ 54f8487c73SHong Zhang #endif 55f8487c73SHong Zhang PetscInt *garray; /* global index of all off-processor columns */ 56f8487c73SHong Zhang 57f8487c73SHong Zhang /* The following variables are used for matrix-vector products */ 58f8487c73SHong Zhang Vec lvec; /* local vector */ 59f8487c73SHong Zhang Vec diag; 60f8487c73SHong Zhang VecScatter Mvctx; /* scatter context for vector */ 61f8487c73SHong Zhang PetscBool roworiented; /* if true, row-oriented input, default true */ 62f8487c73SHong Zhang 63f8487c73SHong Zhang /* The following variables are for MatGetRow() */ 64f8487c73SHong Zhang PetscInt *rowindices; /* column indices for row */ 65f8487c73SHong Zhang PetscScalar *rowvalues; /* nonzero values in row */ 66f8487c73SHong Zhang PetscBool getrowactive; /* indicates MatGetRow(), not restored */ 67f8487c73SHong Zhang 68f8487c73SHong Zhang /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation and nonzero pattern */ 69f8487c73SHong Zhang PetscInt *ld; /* number of entries per row left of diagona block */ 70f8487c73SHong Zhang 71a1a4e84aSHong Zhang /* Used by MatMatMult() and MatPtAP() */ 72f8487c73SHong Zhang Mat_PtAPMPI *ptap; 73ca45077fSPaul Mullowney 74f996eeb8SHong Zhang /* used by MatMatMatMult() */ 75f996eeb8SHong Zhang Mat_MatMatMatMult *matmatmatmult; 76f996eeb8SHong Zhang 77bbf3fe20SPaul Mullowney /* Used by MPICUSP and MPICUSPARSE classes */ 78bbf3fe20SPaul Mullowney void * spptr; 795cc03489SHong Zhang 80f8487c73SHong Zhang } Mat_MPIAIJ; 81f8487c73SHong Zhang 828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MPIAIJ(Mat); 832327d61dSSatish Balay 845a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring); 855a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*); 865a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat); 875a576424SJed Brown PETSC_INTERN PetscErrorCode MatDisAssemble_MPIAIJ(Mat); 885a576424SJed Brown PETSC_INTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat*); 895a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt); 9093dfae19SHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringCreate_MPIXAIJ(Mat,ISColoring,MatFDColoring); 91f86b9fbaSHong Zhang PETSC_INTERN PetscErrorCode MatFDColoringSetUp_MPIXAIJ(Mat,ISColoring,MatFDColoring); 925a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 9353dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSubMatricesMPI_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 945a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat,MatGetSubMatrixOption,MatReuse,Mat *[]); 9553dd7562SDmitry Karpeev 965494a064SHong Zhang 975a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat*); 985a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSubMatrix_MPIAIJ_Private (Mat,IS,IS,PetscInt,MatReuse,Mat*); 993233f0ddSHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,MatReuse,Mat*); 100d6037b41SHong Zhang 1015a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer); 102e3225e9dSHong Zhang PETSC_INTERN PetscErrorCode MatCreateColmap_MPIAIJ_Private(Mat); 1035a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1045a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1050fc8cf34SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 106b2405163SHong Zhang PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1075a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 108f8efd71cSHong Zhang PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 109f8487c73SHong Zhang 1105a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMult_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,MatReuse,PetscReal,Mat*); 1115a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,PetscReal,Mat*); 1125a576424SJed Brown PETSC_INTERN PetscErrorCode MatMatMatMultNumeric_MPIAIJ_MPIAIJ_MPIAIJ(Mat,Mat,Mat,Mat); 113f996eeb8SHong Zhang 1145a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1155a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1165a576424SJed Brown PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 1170d3441aeSHong Zhang 1180d3441aeSHong Zhang PETSC_INTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_new(Mat,Mat,PetscReal,Mat*); 1190d3441aeSHong Zhang PETSC_INTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_new(Mat,Mat,Mat); 1200d3441aeSHong Zhang 1215a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat); 1225a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 123c5af039cSHong Zhang 1245a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 1255a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 1265a576424SJed Brown PETSC_INTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 1275a576424SJed Brown PETSC_INTERN PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*); 1285a576424SJed Brown PETSC_INTERN PetscErrorCode MatSetOption_MPIAIJ(Mat,MatOption,PetscBool); 129187b3c17SHong Zhang 1305a576424SJed Brown PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1312bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,PetscReal,Mat*); 1326da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 1336da69ca6SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 1342bbb1c24SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat,Mat,Mat); 135c216dbf3SHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat,Mat,Mat); 1368949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 1378949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 1388949adfdSHong Zhang PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 1395a576424SJed Brown PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*); 140c0d3702cSSatish Balay 14110e0fddcSSatish Balay PETSC_INTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 1423a7fca6bSBarry Smith 143ce63c4c1SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) 1445a576424SJed Brown PETSC_INTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*); 14597304618SKris Buschelman #endif 1465a576424SJed Brown PETSC_INTERN PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec); 1475a576424SJed Brown PETSC_INTERN PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo*); 14897304618SKris Buschelman 149001ddc4fSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIX_private(PetscInt,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,const PetscInt*,PetscInt*); 150001ddc4fSHong Zhang 15111bd1e4dSLisandro Dalcin extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat*); 1527087cfbeSBarry Smith extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 1533a7fca6bSBarry Smith 15453dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatGetSeqMats_MPIAIJ(Mat,Mat*,Mat*); 15553dd7562SDmitry Karpeev PETSC_INTERN PetscErrorCode MatSetSeqMats_MPIAIJ(Mat,IS,IS,IS,MatStructure,Mat,Mat); 15653dd7562SDmitry Karpeev 157*e497d3c8SHong Zhang /* compute apa = A[i,:]*P = Ad[i,:]*P_loc + Ao*[i,:]*P_oth */ 158*e497d3c8SHong Zhang #define AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa) \ 159*e497d3c8SHong Zhang {\ 160*e497d3c8SHong Zhang PetscInt _anz,_pnz,_j,_k,*_ai,*_aj,_row,*_pi,*_pj; \ 161*e497d3c8SHong Zhang PetscScalar *_aa,_valtmp,*_pa; \ 162*e497d3c8SHong Zhang /* diagonal portion of A */\ 163*e497d3c8SHong Zhang _ai = ad->i;\ 164*e497d3c8SHong Zhang _anz = _ai[i+1] - _ai[i];\ 165*e497d3c8SHong Zhang _aj = ad->j + _ai[i];\ 166*e497d3c8SHong Zhang _aa = ad->a + _ai[i];\ 167*e497d3c8SHong Zhang for (_j=0; _j<_anz; _j++) {\ 168*e497d3c8SHong Zhang _row = _aj[_j]; \ 169*e497d3c8SHong Zhang _pi = p_loc->i; \ 170*e497d3c8SHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 171*e497d3c8SHong Zhang _pj = p_loc->j + _pi[_row]; \ 172*e497d3c8SHong Zhang _pa = p_loc->a + _pi[_row]; \ 173*e497d3c8SHong Zhang /* perform dense axpy */ \ 174*e497d3c8SHong Zhang _valtmp = _aa[_j]; \ 175*e497d3c8SHong Zhang for (_k=0; _k<_pnz; _k++) { \ 176*e497d3c8SHong Zhang apa[_pj[_k]] += _valtmp*_pa[_k]; \ 177*e497d3c8SHong Zhang } \ 178*e497d3c8SHong Zhang PetscLogFlops(2.0*_pnz); \ 179*e497d3c8SHong Zhang } \ 180*e497d3c8SHong Zhang /* off-diagonal portion of A */ \ 181*e497d3c8SHong Zhang _ai = ao->i;\ 182*e497d3c8SHong Zhang _anz = _ai[i+1] - _ai[i]; \ 183*e497d3c8SHong Zhang _aj = ao->j + _ai[i]; \ 184*e497d3c8SHong Zhang _aa = ao->a + _ai[i]; \ 185*e497d3c8SHong Zhang for (_j=0; _j<_anz; _j++) { \ 186*e497d3c8SHong Zhang _row = _aj[_j]; \ 187*e497d3c8SHong Zhang _pi = p_oth->i; \ 188*e497d3c8SHong Zhang _pnz = _pi[_row+1] - _pi[_row]; \ 189*e497d3c8SHong Zhang _pj = p_oth->j + _pi[_row]; \ 190*e497d3c8SHong Zhang _pa = p_oth->a + _pi[_row]; \ 191*e497d3c8SHong Zhang /* perform dense axpy */ \ 192*e497d3c8SHong Zhang _valtmp = _aa[_j]; \ 193*e497d3c8SHong Zhang for (_k=0; _k<_pnz; _k++) { \ 194*e497d3c8SHong Zhang apa[_pj[_k]] += _valtmp*_pa[_k]; \ 195*e497d3c8SHong Zhang } \ 196*e497d3c8SHong Zhang PetscLogFlops(2.0*_pnz); \ 197*e497d3c8SHong Zhang } \ 198*e497d3c8SHong Zhang } 1993369ce9aSBarry Smith #endif 200