1 2 #include <../src/mat/impls/dense/seq/dense.h> 3 4 /* Data stuctures for basic parallel dense matrix */ 5 6 typedef struct { /* used by MatMatMult_MPIDense_MPIDense() */ 7 Mat Ae,Be,Ce; /* matrix in Elemental format */ 8 PetscErrorCode (*destroy)(Mat); 9 } Mat_MatMultDense; 10 11 typedef struct { /* used by MatTransposeMatMult_MPIDense_MPIDense() */ 12 PetscScalar *sendbuf,*atbarray; 13 PetscMPIInt *recvcounts; 14 PetscErrorCode (*destroy)(Mat); 15 PetscMPIInt tag; 16 } Mat_TransMatMultDense; 17 18 typedef struct { /* used by MatMatTransposeMult_MPIDense_MPIDense() */ 19 PetscScalar *buf[2]; 20 PetscMPIInt tag; 21 PetscErrorCode (*destroy)(Mat); 22 } Mat_MatTransMultDense; 23 24 typedef struct { 25 PetscInt nvec; /* this is the n size for the vector one multiplies with */ 26 Mat A; /* local submatrix */ 27 PetscMPIInt size; /* size of communicator */ 28 PetscMPIInt rank; /* rank of proc in communicator */ 29 30 /* The following variables are used for matrix assembly */ 31 PetscBool donotstash; /* Flag indicationg 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 */ 37 38 /* The following variables are used for matrix-vector products */ 39 Vec lvec; /* local vector */ 40 VecScatter Mvctx; /* scatter context for vector */ 41 PetscBool roworiented; /* if true, row oriented input (default) */ 42 43 Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult_MPIAIJ_MPIDense */ 44 Mat_TransMatMultDense *atbdense; /* used by MatTransposeMatMult_MPIDense_MPIDense */ 45 Mat_MatMultDense *abdense; /* used by MatMatMult_MPIDense_MPIDense */ 46 Mat_MatTransMultDense *abtdense; /* used by MatMatTransposeMult_MPIDense_MPIDense */ 47 } Mat_MPIDense; 48 49 PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 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 MatEqual_MPIDense(Mat,Mat,PetscBool*); 53 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 54 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 55 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 56 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 57 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 58 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 59 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 60 61 #if defined(PETSC_HAVE_ELEMENTAL) 62 PETSC_INTERN PetscErrorCode MatMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 63 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 64 #endif 65