Actual source code: mpidense.h

petsc-3.8.4 2018-03-24
Report Typos and Errors

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

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

  6: typedef struct { /* used by MatMatMult_MPIDense_MPIDense() */
  7:   Mat            Ae,Be,Ce;           /* matrix in Elemental format */
  8:   PetscErrorCode (*destroy)(Mat);
  9: } Mat_MatMultDense;

 11: typedef struct { /* used by MatTransposeMatMult_MPIDense_MPIDense() */
 12:   PetscScalar    *sendbuf,*atbarray;
 13:   PetscMPIInt    *recvcounts;
 14:   PetscErrorCode (*destroy)(Mat);
 15: } Mat_TransMatMultDense;

 17: typedef struct {
 18:   PetscInt    nvec;                     /* this is the n size for the vector one multiplies with */
 19:   Mat         A;                        /* local submatrix */
 20:   PetscMPIInt size;                     /* size of communicator */
 21:   PetscMPIInt rank;                     /* rank of proc in communicator */

 23:   /* The following variables are used for matrix assembly */
 24:   PetscBool   donotstash;               /* Flag indicationg if values should be stashed */
 25:   MPI_Request *send_waits;              /* array of send requests */
 26:   MPI_Request *recv_waits;              /* array of receive requests */
 27:   PetscInt    nsends,nrecvs;           /* numbers of sends and receives */
 28:   PetscScalar *svalues,*rvalues;       /* sending and receiving data */
 29:   PetscInt    rmax;                     /* maximum message length */

 31:   /* The following variables are used for matrix-vector products */
 32:   Vec        lvec;                      /* local vector */
 33:   VecScatter Mvctx;                     /* scatter context for vector */
 34:   PetscBool  roworiented;               /* if true, row oriented input (default) */

 36:   Mat_MatTransMatMult   *atb;           /* used by MatTransposeMatMult_MPIAIJ_MPIDense */
 37:   Mat_TransMatMultDense *atbdense;      /* used by MatTransposeMatMult_MPIDense_MPIDense */
 38:   Mat_MatMultDense      *abdense;       /* used by MatMatMult_MPIDense_MPIDense */
 39: } Mat_MPIDense;

 41: PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer);
 42: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
 43: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 44: PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*);
 45: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*);
 46: PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
 47: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
 48: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
 49: PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
 50: PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*);
 51: PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);

 53: #if defined(PETSC_HAVE_ELEMENTAL)
 54: PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
 55: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat);
 56: #endif