Actual source code: mpidense.h

 2:  #include src/mat/impls/dense/seq/dense.h

  4:   /*  Data stuctures for basic parallel dense matrix  */

  6: /* Structure to hold the information for factorization of a dense matrix */
  7: /* Most of this info is used in the pipe send/recv routines */
  8: typedef struct {
  9:   PetscInt    nlnr;        /* number of local rows downstream */
 10:   PetscInt    nrend;       /* rend for downstream processor */
 11:   PetscInt    nbr,pnbr;   /* Down and upstream neighbors */
 12:   PetscInt    *tag;        /* message tags */
 13:   PetscInt    currow;      /* current row number */
 14:   PetscInt    phase;       /* phase (used to indicate tag) */
 15:   PetscInt    up;          /* Are we moving up or down in row number? */
 16:   PetscInt    use_bcast;   /* Are we broadcasting max length? */
 17:   PetscInt    nsend;       /* number of sends */
 18:   PetscInt    nrecv;       /* number of receives */

 20:   /* data initially in matrix context */
 21:   PetscInt    k;           /* Blocking factor (unused as yet) */
 22:   PetscInt    k2;          /* Blocking factor for solves */
 23:   PetscScalar *temp;
 24:   PetscInt    nlptr;
 25:   PetscInt    *lrows;
 26:   PetscInt    *nlrows;
 27:   PetscInt    *pivots;
 28: } FactorCtx;

 30: #define PIPEPHASE (ctx->phase == 0)

 32: typedef struct {
 33:   PetscInt      nvec;                   /* this is the n size for the vector one multiplies with */
 34:   Mat           A;                      /* local submatrix */
 35:   PetscMPIInt   size;                   /* size of communicator */
 36:   PetscMPIInt   rank;                   /* rank of proc in communicator */
 37:   /* The following variables are used for matrix assembly */
 38:   PetscTruth    donotstash;             /* Flag indicationg if values should be stashed */
 39:   MPI_Request   *send_waits;            /* array of send requests */
 40:   MPI_Request   *recv_waits;            /* array of receive requests */
 41:   PetscInt      nsends,nrecvs;         /* numbers of sends and receives */
 42:   PetscScalar   *svalues,*rvalues;     /* sending and receiving data */
 43:   PetscInt      rmax;                   /* maximum message length */

 45:   /* The following variables are used for matrix-vector products */

 47:   Vec           lvec;                   /* local vector */
 48:   VecScatter    Mvctx;                  /* scatter context for vector */

 50:   PetscTruth    roworiented;            /* if true, row oriented input (default) */
 51:   FactorCtx     *factor;                /* factorization context */
 52: } Mat_MPIDense;

 54: EXTERN PetscErrorCode MatLoad_MPIDense(PetscViewer, MatType,Mat*);
 55: EXTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
 56: EXTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 57: EXTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscTruth*);
 58: EXTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*);
 59: EXTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
 60: EXTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
 61: EXTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);