xref: /petsc/src/mat/impls/sell/mpi/mpisell.h (revision 5d7a6ebe9dde080aedbe86be0085708de8f97bb7)
1 #ifndef __MPISELL_H
2 #define __MPISELL_H
3 
4 #include <../src/mat/impls/sell/seq/sell.h>
5 
6 typedef struct {
7   Mat         A, B; /* local submatrices: A (diag part),
8                                            B (off-diag part) */
9   PetscMPIInt size; /* size of communicator */
10   PetscMPIInt rank; /* rank of proc in communicator */
11 
12   /* The following variables are used for matrix assembly */
13   PetscBool    donotstash;        /* PETSC_TRUE if off processor entries dropped */
14   MPI_Request *send_waits;        /* array of send requests */
15   MPI_Request *recv_waits;        /* array of receive requests */
16   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */
17   PetscScalar *svalues, *rvalues; /* sending and receiving data */
18   PetscInt     rmax;              /* maximum message length */
19 #if defined(PETSC_USE_CTABLE)
20   PetscHMapI colmap;
21 #else
22   PetscInt *colmap; /* local col number of off-diag col */
23 #endif
24   PetscInt *garray; /* global index of all off-processor columns */
25 
26   /* The following variables are used for matrix-vector products */
27   Vec        lvec; /* local vector */
28   Vec        diag;
29   VecScatter Mvctx;       /* scatter context for vector */
30   PetscBool  roworiented; /* if true, row-oriented input, default true */
31 
32   /* The following variables are for MatGetRow() */
33   PetscInt    *rowindices;   /* column indices for row */
34   PetscScalar *rowvalues;    /* nonzero values in row */
35   PetscBool    getrowactive; /* indicates MatGetRow(), not restored */
36 
37   PetscInt *ld; /* number of entries per row left of diagona block */
38 } Mat_MPISELL;
39 
40 PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat);
41 
42 PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPISELL(Mat, MatAssemblyType);
43 
44 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPISELL(Mat);
45 PETSC_INTERN PetscErrorCode MatDisAssemble_MPISELL(Mat);
46 
47 PETSC_INTERN PetscErrorCode MatDestroy_MPISELL_PtAP(Mat);
48 PETSC_INTERN PetscErrorCode MatDestroy_MPISELL(Mat);
49 
50 PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPISELL(Mat, Mat *);
51 
52 PETSC_INTERN PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat, MatType, MatReuse, Mat *);
53 PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat, MatType, MatReuse, Mat *);
54 
55 PETSC_INTERN PetscErrorCode MatSOR_MPISELL(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
56 
57 PETSC_INTERN PetscErrorCode MatCreateColmap_MPISELL_Private(Mat);
58 PETSC_INTERN PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat, Vec);
59 #endif
60