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 MatTransposeMatMult_MPIDense_MPIDense() */ 7 PetscScalar *sendbuf,*atbarray; 8 PetscMPIInt *recvcounts; 9 PetscErrorCode (*destroy)(Mat); 10 } Mat_TransMatMultDense; 11 12 typedef struct { 13 PetscInt nvec; /* this is the n size for the vector one multiplies with */ 14 Mat A; /* local submatrix */ 15 PetscMPIInt size; /* size of communicator */ 16 PetscMPIInt rank; /* rank of proc in communicator */ 17 /* The following variables are used for matrix assembly */ 18 PetscBool donotstash; /* Flag indicationg if values should be stashed */ 19 MPI_Request *send_waits; /* array of send requests */ 20 MPI_Request *recv_waits; /* array of receive requests */ 21 PetscInt nsends,nrecvs; /* numbers of sends and receives */ 22 PetscScalar *svalues,*rvalues; /* sending and receiving data */ 23 PetscInt rmax; /* maximum message length */ 24 25 /* The following variables are used for matrix-vector products */ 26 Vec lvec; /* local vector */ 27 VecScatter Mvctx; /* scatter context for vector */ 28 PetscBool roworiented; /* if true, row oriented input (default) */ 29 30 Mat_MatTransMatMult *atb; /* used by MatTransposeMatMult_MPIAIJ_MPIDense */ 31 Mat_TransMatMultDense *atbdense; /* used by MatTransposeMatMult_MPIDense_MPIDense */ 32 } Mat_MPIDense; 33 34 PETSC_INTERN PetscErrorCode MatLoad_MPIDense(Mat,PetscViewer); 35 PETSC_INTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat); 36 PETSC_INTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]); 37 PETSC_INTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscBool*); 38 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 39 PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 40 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*); 41 PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat); 42 PETSC_INTERN PetscErrorCode MatTransposeMatMult_MPIDense_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*); 43 PETSC_INTERN PetscErrorCode MatTransposeMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*); 44 PETSC_INTERN PetscErrorCode MatTransposeMatMultNumeric_MPIDense_MPIDense(Mat,Mat,Mat); 45