xref: /petsc/src/mat/impls/aij/mpi/mpiaij.h (revision b1d57f1590ff1d2f11fead0759f3cc66a5edcec8)
18a729477SBarry Smith 
23369ce9aSBarry Smith #if !defined(__MPIAIJ_H)
33369ce9aSBarry Smith #define __MPIAIJ_H
43369ce9aSBarry Smith 
53369ce9aSBarry Smith #include "src/mat/impls/aij/seq/aij.h"
632b39ab2SSatish Balay #include "src/sys/ctable.h"
78a729477SBarry Smith 
88a729477SBarry Smith typedef struct {
997f1f81fSBarry Smith   PetscInt      *rowners,*cowners;     /* ranges owned by each processor */
1097f1f81fSBarry Smith   PetscInt      rstart,rend;           /* starting and ending owned rows */
1197f1f81fSBarry Smith   PetscInt      cstart,cend;           /* starting and ending owned columns */
128c386bf4SBarry Smith   Mat           A,B;                   /* local submatrices: A (diag part),
138c386bf4SBarry Smith                                            B (off-diag part) */
1497f1f81fSBarry Smith   PetscMPIInt   size;                   /* size of communicator */
1597f1f81fSBarry Smith   PetscMPIInt   rank;                   /* rank of proc in communicator */
168c386bf4SBarry Smith 
178c386bf4SBarry Smith   /* The following variables are used for matrix assembly */
188c386bf4SBarry Smith 
197c922b88SBarry Smith   PetscTruth    donotstash;             /* PETSC_TRUE if off processor entries dropped */
208c386bf4SBarry Smith   MPI_Request   *send_waits;            /* array of send requests */
218c386bf4SBarry Smith   MPI_Request   *recv_waits;            /* array of receive requests */
2297f1f81fSBarry Smith   PetscInt      nsends,nrecvs;         /* numbers of sends and receives */
23ea709b57SSatish Balay   PetscScalar   *svalues,*rvalues;     /* sending and receiving data */
2497f1f81fSBarry Smith   PetscInt      rmax;                   /* maximum message length */
25ac2a4f0dSBarry Smith #if defined (PETSC_USE_CTABLE)
26ac2a4f0dSBarry Smith   PetscTable    colmap;
2732b39ab2SSatish Balay #else
2897f1f81fSBarry Smith   PetscInt      *colmap;                /* local col number of off-diag col */
2932b39ab2SSatish Balay #endif
3097f1f81fSBarry Smith   PetscInt      *garray;                /* global index of all off-processor columns */
318c386bf4SBarry Smith 
328c386bf4SBarry Smith   /* The following variables are used for matrix-vector products */
338c386bf4SBarry Smith 
348c386bf4SBarry Smith   Vec           lvec;              /* local vector */
358c386bf4SBarry Smith   VecScatter    Mvctx;             /* scatter context for vector */
367c922b88SBarry Smith   PetscTruth    roworiented;       /* if true, row-oriented input, default true */
373369ce9aSBarry Smith 
383369ce9aSBarry Smith   /* The following variables are for MatGetRow() */
393369ce9aSBarry Smith 
40*b1d57f15SBarry Smith   PetscInt      *rowindices;       /* column indices for row */
41ea709b57SSatish Balay   PetscScalar   *rowvalues;        /* nonzero values in row */
423369ce9aSBarry Smith   PetscTruth    getrowactive;      /* indicates MatGetRow(), not restored */
43b0a32e0cSBarry Smith 
44180e1086SBarry Smith } Mat_MPIAIJ;
453369ce9aSBarry Smith 
46ff134f7aSHong Zhang typedef struct { /* used by MatMatMult_MPIAIJ_MPIAIJ for reusing symbolic mat product */
47ff134f7aSHong Zhang   IS       isrowa,isrowb,iscolb;
48ff134f7aSHong Zhang   Mat      *aseq,*bseq,C_seq; /* A_seq=aseq[0], B_seq=bseq[0] */
495a5f83e5SHong Zhang   Mat      A_loc,B_seq;
5097f1f81fSBarry Smith   PetscInt brstart; /* starting owned rows of B in matrix bseq[0]; brend = brstart+B->m */
51ff134f7aSHong Zhang } Mat_MatMatMultMPI;
52ff134f7aSHong Zhang 
53b90dcfe3SHong Zhang typedef struct { /* used by MatMerge_SeqsToMPI for reusing the merged matrix */
54b90dcfe3SHong Zhang   PetscMap     rowmap;
55*b1d57f15SBarry Smith   PetscInt     nsend,*bi,*bj,**buf_ri,**buf_rj;
56*b1d57f15SBarry Smith   PetscMPIInt  *len_s,*len_r,*id_r,nrecv; /* array of length of comm->size, store send/recv matrix values */
574087d160SHong Zhang   Mat          C_seq;
58b90dcfe3SHong Zhang } Mat_Merge_SeqsToMPI;
59b90dcfe3SHong Zhang 
60dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetColoring_MPIAIJ(Mat,ISColoring);
61dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValuesAdic_MPIAIJ(Mat,void*);
6297f1f81fSBarry Smith EXTERN PetscErrorCode MatSetValuesAdifor_MPIAIJ(Mat,PetscInt,void*);
63dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetUpMultiply_MPIAIJ(Mat);
64dfbe8321SBarry Smith EXTERN PetscErrorCode DisAssemble_MPIAIJ(Mat);
65dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_MPIAIJ(Mat,MatDuplicateOption,Mat *);
6697f1f81fSBarry Smith EXTERN PetscErrorCode MatIncreaseOverlap_MPIAIJ(Mat,PetscInt,IS [],PetscInt);
67dfbe8321SBarry Smith EXTERN PetscErrorCode MatFDColoringCreate_MPIAIJ(Mat,ISColoring,MatFDColoring);
6897f1f81fSBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_MPIAIJ (Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
6997f1f81fSBarry Smith EXTERN PetscErrorCode MatGetSubMatrix_MPIAIJ (Mat,IS,IS,PetscInt,MatReuse,Mat *);
70dfbe8321SBarry Smith EXTERN PetscErrorCode MatLoad_MPIAIJ(PetscViewer,const MatType,Mat*);
71dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
72dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
73dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
74ff134f7aSHong Zhang EXTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
75ff134f7aSHong Zhang EXTERN PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat,Mat,PetscReal,Mat*);
76ff134f7aSHong Zhang EXTERN PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat,Mat,Mat);
7797f1f81fSBarry Smith EXTERN PetscErrorCode MatSetValues_MPIAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
78c0d3702cSSatish Balay 
793a7fca6bSBarry Smith 
8097304618SKris Buschelman #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
81dfbe8321SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_MPIAIJ_TFS(Mat,IS,IS,MatFactorInfo*,Mat*);
8297304618SKris Buschelman #endif
8397304618SKris Buschelman 
8497304618SKris Buschelman EXTERN_C_BEGIN
85dfbe8321SBarry Smith EXTERN PetscErrorCode MatGetDiagonalBlock_MPIAIJ(Mat,PetscTruth *,MatReuse,Mat *);
86dfbe8321SBarry Smith EXTERN PetscErrorCode MatDiagonalScaleLocal_MPIAIJ(Mat,Vec);
8797304618SKris Buschelman EXTERN_C_END
883a7fca6bSBarry Smith 
893369ce9aSBarry Smith #endif
90