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