Actual source code: matisimpl.h

  1: #pragma once

  3: #include <petscsf.h>
  4: #include <petsc/private/matimpl.h>

  6: typedef struct {
  7:   Mat          A;       /* the local matrix */
  8:   VecScatter   cctx;    /* column scatter */
  9:   VecScatter   rctx;    /* row scatter */
 10:   Vec          x, y;    /* work space for matrix vector product */
 11:   Vec          counter; /* counter vector */
 12:   PetscBool    pure_neumann;
 13:   PetscSF      sf, csf; /* SFs for rows and cols */
 14:   PetscInt    *sf_rootdata, *sf_leafdata;
 15:   PetscInt    *csf_rootdata, *csf_leafdata;
 16:   IS           getsub_ris, getsub_cis; /* row and column ISs for MatCreateSubMatrix and MAT_REUSE_MATRIX */
 17:   PetscBool    allow_repeated;         /* allow repeated entries in the local to global map */
 18:   PetscBool    islocalref;             /* is a reference to a local submatrix? */
 19:   PetscBool    locempty;               /* adapt local matrices for empty rows/cols during MatAssemblyEnd_IS */
 20:   PetscBool    storel2l;               /* carry over local-to-local inherited in MatPtAP */
 21:   char        *lmattype;
 22:   PetscScalar *bdiag; /* Used by MatInvertBlockDiagonal_IS */

 24:   PetscObjectState lnnzstate; /* nonzero state of local matrix */

 26:   PetscBool keepassembled; /* store assembled form if needed */
 27:   Mat       assembledA;    /* assembled operator */
 28:   Mat       dA;            /* For MatGetDiagonalBlock_IS */

 30:   /* Support for negative or repeated entries in l2map
 31:      These maps can be different than the ones passed in by the user via
 32:      MatSetLocalToGlobalMapping */
 33:   ISLocalToGlobalMapping rmapping, cmapping;
 34: } Mat_IS;

 36: struct _MatISLocalFields {
 37:   PetscInt nr, nc;
 38:   IS      *rf, *cf;
 39: };
 40: typedef struct _MatISLocalFields *MatISLocalFields;

 42: struct _MatISPtAP {
 43:   PetscReal fill;
 44:   IS        cis0, cis1, ris0, ris1;
 45:   Mat      *lP;
 46: };
 47: typedef struct _MatISPtAP *MatISPtAP;

 49: PETSC_EXTERN PetscErrorCode MatISSetMPIXAIJPreallocation_Private(Mat, Mat, PetscBool);