1 2 #if !defined(__MPIAIJ_H) 3 #define __MPIAIJ_H 4 5 #include "src/mat/impls/aij/seq/aij.h" 6 #include "src/sys/ctable.h" 7 8 typedef struct { 9 Mat A,B; /* local submatrices: A (diag part), 10 B (off-diag part) */ 11 PetscMPIInt size; /* size of communicator */ 12 PetscMPIInt rank; /* rank of proc in communicator */ 13 14 /* The following variables are used for matrix assembly */ 15 16 PetscTruth donotstash; /* PETSC_TRUE if off processor entries dropped */ 17 MPI_Request *send_waits; /* array of send requests */ 18 MPI_Request *recv_waits; /* array of receive requests */ 19 PetscInt nsends,nrecvs; /* numbers of sends and receives */ 20 PetscScalar *svalues,*rvalues; /* sending and receiving data */ 21 PetscInt rmax; /* maximum message length */ 22 #if defined (PETSC_USE_CTABLE) 23 PetscTable colmap; 24 #else 25 PetscInt *colmap; /* local col number of off-diag col */ 26 #endif 27 PetscInt *garray; /* global index of all off-processor columns */ 28 29 /* The following variables are used for matrix-vector products */ 30 31 Vec lvec; /* local vector */ 32 VecScatter Mvctx; /* scatter context for vector */ 33 PetscTruth roworiented; /* if true, row-oriented input, default true */ 34 35 /* The following variables are for MatGetRow() */ 36 37 PetscInt *rowindices; /* column indices for row */ 38 PetscScalar *rowvalues; /* nonzero values in row */ 39 PetscTruth getrowactive; /* indicates MatGetRow(), not restored */ 40 41 } Mat_MPIAIJ; 42 43 typedef struct { /* used by MatMatMult_MPIAIJ_MPIAIJ and MatPtAP_MPIAIJ_MPIAIJ for reusing symbolic mat product */ 44 PetscInt *startsj; 45 PetscScalar *bufa; 46 IS isrowa,isrowb,iscolb; 47 Mat *aseq,*bseq,C_seq; /* A_seq=aseq[0], B_seq=bseq[0] */ 48 Mat A_loc,B_seq; 49 Mat B_loc,B_oth; /* partial B_seq -- intend to replace B_seq */ 50 PetscInt brstart; /* starting owned rows of B in matrix bseq[0]; brend = brstart+B->m */ 51 PetscInt *abi,*abj; /* symbolic i and j arrays of the local product A_loc*B_seq */ 52 PetscInt abnz_max; /* max(abi[i+1] - abi[i]), max num of nnz in a row of A_loc*B_seq */ 53 MatReuse reuse; 54 PetscErrorCode (*MatDestroy)(Mat); 55 } Mat_MatMatMultMPI; 56 57 typedef struct { /* used by MatMerge_SeqsToMPI for reusing the merged matrix */ 58 PetscMap rowmap; 59 PetscInt **buf_ri,**buf_rj; 60 PetscMPIInt *len_s,*len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 61 PetscMPIInt nsend,nrecv; 62 PetscInt *bi,*bj; /* i and j array of the local portion of mpi C=P^T*A*P */ 63 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 64 } Mat_Merge_SeqsToMPI; 65 66 EXTERN PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring); 67 EXTERN PetscErrorCode MatSetValuesAdic_MPIAIJ(Mat,void*); 68 EXTERN PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*); 69 EXTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat); 70 EXTERN PetscErrorCode DisAssemble_MPIAIJ(Mat); 71 EXTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat *); 72 EXTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt); 73 EXTERN PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 74 EXTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 75 EXTERN PetscErrorCode MatGetSubMatrix_MPIAIJ (Mat,IS,IS,PetscInt,MatReuse,Mat *); 76 EXTERN PetscErrorCode MatLoad_MPIAIJ(PetscViewer, MatType,Mat*); 77 EXTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 78 EXTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 79 EXTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 80 EXTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat,Mat,PetscReal,Mat*); 81 EXTERN PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat,Mat,Mat); 82 EXTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 83 EXTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 84 EXTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 85 EXTERN PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 86 EXTERN PetscErrorCode PetscObjectContainerDestroy_Mat_MatMatMultMPI(void*); 87 88 EXTERN_C_BEGIN 89 EXTERN PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 90 EXTERN_C_END 91 92 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) 93 EXTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,MatFactorInfo*,Mat*); 94 #endif 95 96 EXTERN_C_BEGIN 97 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatGetDiagonalBlock_MPIAIJ(Mat,PetscTruth *,MatReuse,Mat *); 98 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 99 EXTERN_C_END 100 101 #endif 102