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