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 MatMatMult_MPIAIJ_MPIAIJ_32 - implementation used in PETSc-3.2 */ 8 /* used by MatMatMult_MPIAIJ_MPIAIJ() */ 9 IS isrowa,isrowb,iscolb; 10 Mat B_seq,A_loc,C_seq; 11 PetscErrorCode (*destroy)(Mat); 12 PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 13 } Mat_MatMatMultMPI; 14 15 typedef struct { /* used by MatMerge_SeqsToMPI for reusing the merged matrix - implementation used in PETSc-3.2 */ 16 PetscLayout rowmap; 17 PetscInt **buf_ri,**buf_rj; 18 PetscMPIInt *len_s,*len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 19 PetscMPIInt nsend,nrecv; 20 PetscInt *bi,*bj; /* i and j array of the local portion of mpi C=P^T*A*P */ 21 PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 22 PetscErrorCode (*destroy)(Mat); 23 PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 24 } Mat_Merge_SeqsToMPI; 25 26 typedef struct { /* used by MatPtAP_MPIAIJ_MPIAIJ() and MatMatMult_MPIAIJ_MPIAIJ() */ 27 PetscInt *startsj,*startsj_r; 28 PetscScalar *bufa; 29 Mat P_loc,P_oth; /* partial B_seq -- intend to replace B_seq */ 30 PetscInt *api,*apj; /* symbolic i and j arrays of the local product A_loc*B_seq */ 31 PetscInt abnz_max; /* max(abi[i+1] - abi[i]), max num of nnz in a row of A_loc*B_seq */ 32 MatReuse reuse; 33 PetscScalar *apa; /* tmp array for store a row of A*P used in MatMatMult() */ 34 35 Mat_Merge_SeqsToMPI *merge; 36 PetscErrorCode (*destroy)(Mat); 37 PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*); 38 } Mat_PtAPMPI; 39 40 typedef struct { 41 Mat A,B; /* local submatrices: A (diag part), 42 B (off-diag part) */ 43 PetscMPIInt size; /* size of communicator */ 44 PetscMPIInt rank; /* rank of proc in communicator */ 45 46 /* The following variables are used for matrix assembly */ 47 PetscBool donotstash; /* PETSC_TRUE if off processor entries dropped */ 48 MPI_Request *send_waits; /* array of send requests */ 49 MPI_Request *recv_waits; /* array of receive requests */ 50 PetscInt nsends,nrecvs; /* numbers of sends and receives */ 51 PetscScalar *svalues,*rvalues; /* sending and receiving data */ 52 PetscInt rmax; /* maximum message length */ 53 #if defined (PETSC_USE_CTABLE) 54 PetscTable colmap; 55 #else 56 PetscInt *colmap; /* local col number of off-diag col */ 57 #endif 58 PetscInt *garray; /* global index of all off-processor columns */ 59 60 /* The following variables are used for matrix-vector products */ 61 Vec lvec; /* local vector */ 62 Vec diag; 63 VecScatter Mvctx; /* scatter context for vector */ 64 PetscBool roworiented; /* if true, row-oriented input, default true */ 65 66 /* The following variables are for MatGetRow() */ 67 PetscInt *rowindices; /* column indices for row */ 68 PetscScalar *rowvalues; /* nonzero values in row */ 69 PetscBool getrowactive; /* indicates MatGetRow(), not restored */ 70 71 /* Used by MatDistribute_MPIAIJ() to allow reuse of previous matrix allocation and nonzero pattern */ 72 PetscInt *ld; /* number of entries per row left of diagona block */ 73 74 /* Used by MatMatMult() and MatPtAP() */ 75 Mat_PtAPMPI *ptap; 76 } Mat_MPIAIJ; 77 78 extern PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring); 79 extern PetscErrorCode MatSetValuesAdic_MPIAIJ(Mat,void*); 80 extern PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*); 81 extern PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat); 82 extern PetscErrorCode DisAssemble_MPIAIJ(Mat); 83 extern PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat *); 84 extern PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt); 85 extern PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring); 86 extern PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 87 extern PetscErrorCode MatGetSubMatrix_MPIAIJ_All(Mat,MatGetSubMatrixOption,MatReuse,Mat *[]); 88 extern PetscErrorCode MatGetSubMatricesParallel_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 89 90 extern PetscErrorCode MatGetSubMatrix_MPIAIJ(Mat,IS,IS,MatReuse,Mat *); 91 extern PetscErrorCode MatGetSubMatrix_MPIAIJ_Private (Mat,IS,IS,PetscInt,MatReuse,Mat *); 92 extern PetscErrorCode MatGetMultiProcBlock_MPIAIJ(Mat,MPI_Comm,Mat*); 93 94 extern PetscErrorCode MatLoad_MPIAIJ(Mat,PetscViewer); 95 extern PetscErrorCode MatMatMult_MPIDense_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 96 extern PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 97 extern PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 98 extern PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_32(Mat,Mat,PetscReal,Mat*); 99 extern PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 100 101 extern PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat,Mat,PetscReal,Mat*); 102 extern PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat,Mat,Mat); 103 extern PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*); 104 extern PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat); 105 extern PetscErrorCode MatGetBrowsOfAoCols_MPIAIJ(Mat,Mat,MatReuse,PetscInt**,PetscInt**,MatScalar**,Mat*); 106 extern PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode); 107 extern PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat); 108 extern PetscErrorCode MatDestroy_MPIAIJ_MatMatMult_32(Mat); 109 extern PetscErrorCode PetscContainerDestroy_Mat_MatMatMultMPI(void*); 110 extern PetscErrorCode MatGetRedundantMatrix_MPIAIJ(Mat,PetscInt,MPI_Comm,PetscInt,MatReuse,Mat*); 111 extern PetscErrorCode MatGetSeqNonzeroStructure_MPIAIJ(Mat,Mat*); 112 113 EXTERN_C_BEGIN 114 extern PetscErrorCode MatMPIAIJSetPreallocation_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]); 115 EXTERN_C_END 116 117 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_REAL_SINGLE) && !defined(PETSC_USE_REAL___FLOAT128) 118 extern PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,const MatFactorInfo*,Mat*); 119 #endif 120 extern PetscErrorCode MatSolve_MPIAIJ(Mat,Vec,Vec); 121 extern PetscErrorCode MatILUFactor_MPIAIJ(Mat,IS,IS,const MatFactorInfo *); 122 123 EXTERN_C_BEGIN 124 extern PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,Mat *); 125 extern PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec); 126 EXTERN_C_END 127 128 #endif 129