Actual source code: mpisell.h

  1: #pragma once

  3: #include <../src/mat/impls/sell/seq/sell.h>

  5: typedef struct {
  6:   Mat         A, B; /* local submatrices: A (diag part),
  7:                                            B (off-diag part) */
  8:   PetscMPIInt size; /* size of communicator */
  9:   PetscMPIInt rank; /* rank of proc in communicator */

 11:   /* The following variables are used for matrix assembly */
 12:   PetscBool    donotstash;        /* PETSC_TRUE if off processor entries dropped */
 13:   MPI_Request *send_waits;        /* array of send requests */
 14:   MPI_Request *recv_waits;        /* array of receive requests */
 15:   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */
 16:   PetscScalar *svalues, *rvalues; /* sending and receiving data */
 17:   PetscInt     rmax;              /* maximum message length */
 18: #if defined(PETSC_USE_CTABLE)
 19:   PetscHMapI colmap;
 20: #else
 21:   PetscInt *colmap; /* local col number of off-diag col */
 22: #endif
 23:   PetscInt *garray; /* global index of all off-processor columns */

 25:   /* The following variables are used for matrix-vector products */
 26:   Vec        lvec; /* local vector */
 27:   Vec        diag;
 28:   VecScatter Mvctx;       /* scatter context for vector */
 29:   PetscBool  roworiented; /* if true, row-oriented input, default true */

 31:   /* The following variables are for MatGetRow() */
 32:   PetscInt    *rowindices;   /* column indices for row */
 33:   PetscScalar *rowvalues;    /* nonzero values in row */
 34:   PetscBool    getrowactive; /* indicates MatGetRow(), not restored */

 36:   PetscInt *ld; /* number of entries per row left of diagona block */
 37: } Mat_MPISELL;

 39: PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat);

 41: PETSC_INTERN PetscErrorCode MatAssemblyEnd_MPISELL(Mat, MatAssemblyType);

 43: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPISELL(Mat);
 44: PETSC_INTERN PetscErrorCode MatDisAssemble_MPISELL(Mat);

 46: PETSC_INTERN PetscErrorCode MatDestroy_MPISELL_PtAP(Mat);
 47: PETSC_INTERN PetscErrorCode MatDestroy_MPISELL(Mat);

 49: PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPISELL(Mat, Mat *);

 51: PETSC_INTERN PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat, MatType, MatReuse, Mat *);
 52: PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat, MatType, MatReuse, Mat *);

 54: PETSC_INTERN PetscErrorCode MatSOR_MPISELL(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);

 56: PETSC_INTERN PetscErrorCode MatCreateColmap_MPISELL_Private(Mat);
 57: PETSC_INTERN PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat, Vec);