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