xref: /petsc/src/mat/impls/dense/mpi/mpidense.h (revision 01b828861bdb7b25615807b947b5d42397386d0c)
18965ea79SLois Curfman McInnes 
270f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h"
38965ea79SLois Curfman McInnes 
4052cd425SLois Curfman McInnes   /*  Data stuctures for basic parallel dense matrix  */
5052cd425SLois Curfman McInnes 
68965ea79SLois Curfman McInnes typedef struct {
713f74950SBarry Smith   PetscInt      nvec;                   /* this is the n size for the vector one multiplies with */
8bad5519cSBarry Smith   Mat           A;                      /* local submatrix */
913f74950SBarry Smith   PetscMPIInt   size;                   /* size of communicator */
1013f74950SBarry Smith   PetscMPIInt   rank;                   /* rank of proc in communicator */
118965ea79SLois Curfman McInnes   /* The following variables are used for matrix assembly */
12273d9f13SBarry Smith   PetscTruth    donotstash;             /* Flag indicationg if values should be stashed */
138965ea79SLois Curfman McInnes   MPI_Request   *send_waits;            /* array of send requests */
148965ea79SLois Curfman McInnes   MPI_Request   *recv_waits;            /* array of receive requests */
1513f74950SBarry Smith   PetscInt      nsends,nrecvs;         /* numbers of sends and receives */
16ea709b57SSatish Balay   PetscScalar   *svalues,*rvalues;     /* sending and receiving data */
1713f74950SBarry Smith   PetscInt      rmax;                   /* maximum message length */
188965ea79SLois Curfman McInnes 
198965ea79SLois Curfman McInnes   /* The following variables are used for matrix-vector products */
208965ea79SLois Curfman McInnes 
218965ea79SLois Curfman McInnes   Vec           lvec;                   /* local vector */
22bad5519cSBarry Smith   VecScatter    Mvctx;                  /* scatter context for vector */
23052cd425SLois Curfman McInnes 
24273d9f13SBarry Smith   PetscTruth    roworiented;            /* if true, row oriented input (default) */
258965ea79SLois Curfman McInnes } Mat_MPIDense;
2697304618SKris Buschelman 
27f69a0ea3SMatthew Knepley EXTERN PetscErrorCode MatLoad_MPIDense(PetscViewer, MatType,Mat*);
28dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetUpMultiply_MPIDense(Mat);
2913f74950SBarry Smith EXTERN PetscErrorCode MatGetSubMatrices_MPIDense(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
306e4ee0c6SHong Zhang EXTERN PetscErrorCode MatEqual_MPIDense(Mat,Mat,PetscTruth*);
314ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMultSymbolic_MPIDense_MPIDense(Mat,Mat,PetscReal,Mat*);
324ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat,Mat,MatReuse,PetscReal,Mat*);
334ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat*);
344ae313f4SHong Zhang EXTERN PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
35*01b82886SBarry Smith 
36*01b82886SBarry Smith #if defined(PETSC_HAVE_PLAPACK)
37*01b82886SBarry Smith EXTERN_C_BEGIN
38*01b82886SBarry Smith #include "PLA.h"
39*01b82886SBarry Smith EXTERN_C_END
40*01b82886SBarry Smith 
41*01b82886SBarry Smith typedef struct {
42*01b82886SBarry Smith   MPI_Comm       comm_2d;
43*01b82886SBarry Smith   PLA_Obj        A,pivots;
44*01b82886SBarry Smith   PLA_Template   templ;
45*01b82886SBarry Smith   MPI_Datatype   datatype;
46*01b82886SBarry Smith   PetscInt       nb,nb_alg,ierror,rstart;
47*01b82886SBarry Smith   VecScatter     ctx;
48*01b82886SBarry Smith   IS             is_pla,is_petsc;
49*01b82886SBarry Smith   PetscTruth     pla_solved;
50*01b82886SBarry Smith   MatStructure   mstruct;
51*01b82886SBarry Smith   PetscMPIInt    nprows,npcols;
52*01b82886SBarry Smith 
53*01b82886SBarry Smith   /* Flag to clean up (non-global) Plapack objects during Destroy */
54*01b82886SBarry Smith   PetscTruth CleanUpPlapack;
55*01b82886SBarry Smith } Mat_Plapack;
56*01b82886SBarry Smith 
57*01b82886SBarry Smith #endif
58