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 MatMatMultxxx_MPIDense_MPIDense() */ 7 Mat Ae,Be,Ce; /* matrix in Elemental format */ 8 PetscErrorCode (*destroy)(Mat); 9 } Mat_MatMultDense; 10 11 typedef struct { /* used by MatTransposeMatMultXXX_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 MatMatTransposeMultxxx_MPIDense_MPIDense() */ 19 PetscScalar *buf[2]; 20 PetscMPIInt tag; 21 PetscMPIInt *recvcounts; 22 PetscMPIInt *recvdispls; 23 PetscErrorCode (*destroy)(Mat); 24 PetscInt alg; /* algorithm used */ 25 } Mat_MatTransMultDense; 26 27 typedef struct { 28 PetscInt nvec; /* this is the n size for the vector one multiplies with */ 29 Mat A; /* local submatrix */ 30 PetscMPIInt size; /* size of communicator */ 31 PetscMPIInt rank; /* rank of proc in communicator */ 32 33 /* The following variables are used for matrix assembly */ 34 PetscBool donotstash; /* Flag indicationg if values should be stashed */ 35 MPI_Request *send_waits; /* array of send requests */ 36 MPI_Request *recv_waits; /* array of receive requests */ 37 PetscInt nsends,nrecvs; /* numbers of sends and receives */ 38 PetscScalar *svalues,*rvalues; /* sending and receiving data */ 39 PetscInt rmax; /* maximum message length */ 40 41 /* The following variables are used for matrix-vector products */ 42 Vec lvec; /* local vector */ 43 VecScatter Mvctx; /* scatter context for vector */ 44 PetscBool roworiented; /* if true, row oriented input (default) */ 45 46 Mat_MatTransMatMult *atb; /* used by MatTransposeMatMultxxx_MPIAIJ_MPIDense */ 47 Mat_TransMatMultDense *atbdense; /* used by MatTransposeMatMultxxx_MPIDense_MPIDense */ 48 Mat_MatMultDense *abdense; /* used by MatMatMultxxx_MPIDense_MPIDense */ 49 Mat_MatTransMultDense *abtdense; /* used by MatMatTransposeMultxxx_MPIDense_MPIDense */ 50 } Mat_MPIDense; 51 52 PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 53 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 54 PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 55 PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*); 56 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat); 57 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat); 58 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat); 59 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 60 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat); 61 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 62 63 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense_MPIAIJ(Mat); 64 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIDense(Mat); 65 66 #if defined(PETSC_HAVE_ELEMENTAL) 67 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 68 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Elemental(Mat,Mat,PetscReal,Mat); 69 PETSC_INTERN PetscErrorCode MatMatMultNumeric_Elemental(Mat,Mat,Mat); 70 #endif 71