Actual source code: mpidense.h


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

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

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

 11: typedef struct { /* used by MatTransposeMatMultXXX_MPIDense_MPIDense() */
 12:   PetscScalar *sendbuf;
 13:   Mat         atb;
 14:   PetscMPIInt *recvcounts;
 15:   PetscMPIInt tag;
 16: } Mat_TransMatMultDense;

 18: typedef struct { /* used by MatMatTransposeMultxxx_MPIDense_MPIDense() */
 19:   PetscScalar *buf[2];
 20:   PetscMPIInt tag;
 21:   PetscMPIInt *recvcounts;
 22:   PetscMPIInt *recvdispls;
 23:   PetscInt    alg; /* algorithm used */
 24: } Mat_MatTransMultDense;

 26: typedef struct {
 27:   Mat         A;                        /* local submatrix */

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

 37:   /* The following variables are used for matrix-vector products */
 38:   Vec        lvec;                      /* local vector */
 39:   PetscSF    Mvctx;                     /* for mat-mult communications */
 40:   PetscBool  roworiented;               /* if true, row oriented input (default) */

 42:   /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
 43:   Mat               cmat;      /* matrix representation of a given subset of columns */
 44:   Vec               cvec;      /* vector representation of a given column */
 45:   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
 46:   PetscInt          vecinuse;  /* if cvec is in use (col = vecinuse-1) */
 47:   PetscInt          matinuse;  /* if cmat is in use (cbegin = matinuse-1) */
 48: } Mat_MPIDense;

 50: PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
 51: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 52: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat);

 54: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat);
 55: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat);

 57: #if defined(PETSC_HAVE_ELEMENTAL)
 58: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat,Mat,PetscReal,Mat);
 59: PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat,Mat,Mat);
 60: #endif
 61: PETSC_INTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat);