Actual source code: mpisell.h

petsc-3.14.6 2021-03-30
Report Typos and Errors

  4: #endif
  5: #include <../src/mat/impls/sell/seq/sell.h>

  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 */

 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:   PetscTable 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 */

 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 */

 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 */

 38:   /* Used by MatDistribute_MPISELL() to allow reuse of previous matrix allocation  and nonzero pattern */
 39:   PetscInt *ld;                    /* number of entries per row left of diagona block */
 40: } Mat_MPISELL;

 42: PETSC_EXTERN PetscErrorCode MatCreate_MPISELL(Mat);
 43: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPISELL(Mat);

 45: PETSC_INTERN PetscErrorCode MatDisAssemble_MPISELL(Mat);
 46: PETSC_INTERN PetscErrorCode MatDuplicate_MPISELL(Mat,MatDuplicateOption,Mat*);

 48: PETSC_INTERN PetscErrorCode MatDestroy_MPISELL_PtAP(Mat);
 49: PETSC_INTERN PetscErrorCode MatDestroy_MPISELL(Mat);

 51: PETSC_INTERN PetscErrorCode MatSetValues_MPISELL(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar [],InsertMode);
 52: PETSC_INTERN PetscErrorCode MatSetOption_MPISELL(Mat,MatOption,PetscBool);
 53: PETSC_INTERN PetscErrorCode MatGetSeqNonzeroStructure_MPISELL(Mat,Mat*);

 55: PETSC_INTERN PetscErrorCode MatSetFromOptions_MPISELL(PetscOptionItems*,Mat);
 56: PETSC_INTERN PetscErrorCode MatMPISELLSetPreallocation_MPISELL(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]);

 58: PETSC_INTERN PetscErrorCode MatConvert_MPISELL_MPIAIJ(Mat,MatType,MatReuse,Mat*);
 59: PETSC_INTERN PetscErrorCode MatConvert_MPIAIJ_MPISELL(Mat,MatType,MatReuse,Mat*);

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

 63: PETSC_INTERN PetscErrorCode MatCreateColmap_MPISELL_Private(Mat);
 64: PETSC_INTERN PetscErrorCode MatDiagonalScaleLocal_MPISELL(Mat,Vec);