Actual source code: mpidense.h

  1: #pragma once

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

  6: /*  Data structures for basic parallel dense matrix  */

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

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

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

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

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

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

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

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

 55: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(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);

 63: PETSC_INTERN PetscErrorCode MatShift_MPIDense(Mat, PetscScalar);
 64: PETSC_INTERN PetscErrorCode MatDenseGetColumnVecWrite_MPIDense(Mat, PetscInt, Vec *);
 65: PETSC_INTERN PetscErrorCode MatDenseGetColumnVecRead_MPIDense(Mat, PetscInt, Vec *);
 66: PETSC_INTERN PetscErrorCode MatDenseGetColumnVec_MPIDense(Mat, PetscInt, Vec *);
 67: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecWrite_MPIDense(Mat, PetscInt, Vec *);
 68: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVecRead_MPIDense(Mat, PetscInt, Vec *);
 69: PETSC_INTERN PetscErrorCode MatDenseRestoreColumnVec_MPIDense(Mat, PetscInt, Vec *);

 71: PETSC_INTERN PetscErrorCode MatCreate_MPIDense(Mat);
 72: PETSC_INTERN PetscErrorCode MatGetDiagonal_MPIDense(Mat, Vec);

 74: #if PetscDefined(HAVE_CUDA)
 75: PETSC_INTERN PetscErrorCode MatConvert_MPIDense_MPIDenseCUDA(Mat, MatType, MatReuse, Mat *);
 76: #endif

 78: #if PetscDefined(HAVE_HIP)
 79: PETSC_INTERN PetscErrorCode MatConvert_MPIDense_MPIDenseHIP(Mat, MatType, MatReuse, Mat *);
 80: #endif