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