xref: /petsc/src/mat/impls/baij/mpi/mpibaij.h (revision fdfbdca60871d1230bda7d5d096e2bd2d2a23f7a)
1a4963045SJacob Faibussowitsch #pragma once
2c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h>
3001ddc4fSHong Zhang #include <../src/mat/impls/aij/mpi/mpiaij.h>
4eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h>
579bdfe76SSatish Balay 
6aa482453SBarry Smith #if defined(PETSC_USE_CTABLE)
7eec179cfSJacob Faibussowitsch   #define PETSCTABLE PetscHMapI
832b39ab2SSatish Balay #else
9b24ad042SBarry Smith   #define PETSCTABLE PetscInt *
1032b39ab2SSatish Balay #endif
11e8823a70SSatish Balay 
126c6c5352SBarry Smith #define MPIBAIJHEADER \
13d0f46423SBarry Smith   PetscInt   *rangebs;                            /* rmap->range/bs */ \
14899cda47SBarry Smith   PetscInt    rstartbs, rendbs, cstartbs, cendbs; /* map values / bs  */ \
151d18e487SKris Buschelman   Mat         A, B;                               /* local submatrices: A (diag part), B (off-diag part) */ \
16b24ad042SBarry Smith   PetscMPIInt size;                               /* size of communicator */ \
17b24ad042SBarry Smith   PetscMPIInt rank;                               /* rank of proc in communicator */ \
18521d7252SBarry Smith   PetscInt    bs2;                                /* block size, bs2 = bs*bs */ \
19b24ad042SBarry Smith   PetscInt    Mbs, Nbs;                           /* number block rows/cols in matrix; M/bs, N/bs */ \
20b24ad042SBarry Smith   PetscInt    mbs, nbs;                           /* number block rows/cols on processor; m/bs, n/bs */ \
216c6c5352SBarry Smith \
226c6c5352SBarry Smith   /* The following variables are used for matrix assembly */ \
236c6c5352SBarry Smith \
24ace3abfcSBarry Smith   PetscBool    donotstash;              /* if 1, off processor entries dropped */ \
2597da8949SJed Brown   PetscBool    subset_off_proc_entries; /* PETSC_TRUE if assembly will always communicate a subset of the entries communicated the first time */ \
266c6c5352SBarry Smith   MPI_Request *send_waits;              /* array of send requests */ \
276c6c5352SBarry Smith   MPI_Request *recv_waits;              /* array of receive requests */ \
28b24ad042SBarry Smith   PetscInt     nsends, nrecvs;          /* numbers of sends and receives */ \
296c6c5352SBarry Smith   MatScalar   *svalues, *rvalues;       /* sending and receiving data */ \
30b24ad042SBarry Smith   PetscInt     rmax;                    /* maximum message length */ \
316c6c5352SBarry Smith   PETSCTABLE   colmap;                  /* local col number of off-diag col */ \
326c6c5352SBarry Smith \
33b24ad042SBarry Smith   PetscInt *garray; /* work array */ \
346c6c5352SBarry Smith \
356c6c5352SBarry Smith   /* The following variable is used by blocked matrix assembly */ \
366c6c5352SBarry Smith   MatScalar *barray; /* Block array of size bs2 */ \
376c6c5352SBarry Smith \
386c6c5352SBarry Smith   /* The following variables are used for matrix-vector products */ \
396c6c5352SBarry Smith \
406c6c5352SBarry Smith   Vec        lvec;        /* local vector */ \
416c6c5352SBarry Smith   VecScatter Mvctx;       /* scatter context for vector */ \
42ace3abfcSBarry Smith   PetscBool  roworiented; /* if true, row-oriented input, default true */ \
436c6c5352SBarry Smith \
446c6c5352SBarry Smith   /* The following variables are for MatGetRow() */ \
456c6c5352SBarry Smith \
46b24ad042SBarry Smith   PetscInt    *rowindices;   /* column indices for row */ \
476c6c5352SBarry Smith   PetscScalar *rowvalues;    /* nonzero values in row */ \
48ace3abfcSBarry Smith   PetscBool    getrowactive; /* indicates MatGetRow(), not restored */ \
496c6c5352SBarry Smith \
506c6c5352SBarry Smith   /* Some variables to make MatSetValues and others more efficient */ \
51b24ad042SBarry Smith   PetscInt    rstart_bs, rend_bs; \
52b24ad042SBarry Smith   PetscInt    cstart_bs, cend_bs; \
53b24ad042SBarry Smith   PetscInt   *ht; /* Hash table to speed up matrix assembly */ \
546c6c5352SBarry Smith   MatScalar **hd; /* Hash table data */ \
55b24ad042SBarry Smith   PetscInt    ht_size; \
56b24ad042SBarry Smith   PetscInt    ht_total_ct, ht_insert_ct; /* Hash table statistics */ \
57ace3abfcSBarry Smith   PetscBool   ht_flag;                   /* Flag to indicate if hash tables are used */ \
586c6c5352SBarry Smith   double      ht_fact;                   /* Factor to determine the HT size */ \
596c6c5352SBarry Smith \
60b24ad042SBarry Smith   PetscInt   setvalueslen;  /* only used for single precision computations */ \
617a868f3eSHong Zhang   MatScalar *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied*/ \
627a868f3eSHong Zhang                             /* before calling MatSetValuesXXX_MPIBAIJ_MatScalar() */ \
6326cec326SBarry Smith   PetscBool      ijonly;    /* used in  MatCreateSubMatrices_MPIBAIJ_local() for getting ij structure only */ \
6426cec326SBarry Smith   struct _MatOps cops
65ab9863d7SBarry Smith 
666c6c5352SBarry Smith typedef struct {
671d18e487SKris Buschelman   MPIBAIJHEADER;
68e8823a70SSatish Balay } Mat_MPIBAIJ;
6979bdfe76SSatish Balay 
70b51a4376SLisandro Dalcin PETSC_INTERN PetscErrorCode MatView_MPIBAIJ(Mat, PetscViewer);
715a576424SJed Brown PETSC_INTERN PetscErrorCode MatLoad_MPIBAIJ(Mat, PetscViewer);
72b51a4376SLisandro Dalcin PETSC_INTERN PetscErrorCode MatView_MPIBAIJ_Binary(Mat, PetscViewer);
73b51a4376SLisandro Dalcin PETSC_INTERN PetscErrorCode MatLoad_MPIBAIJ_Binary(Mat, PetscViewer);
74b51a4376SLisandro Dalcin 
755a576424SJed Brown PETSC_INTERN PetscErrorCode MatCreateColmap_MPIBAIJ_Private(Mat);
767dae84e0SHong Zhang PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
77*fdfbdca6SPierre Jolivet PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *, PetscBool);
78*fdfbdca6SPierre Jolivet PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIBAIJ_Private(Mat, IS, IS, PetscInt, MatReuse, Mat *, PetscBool);
79e8271787SHong Zhang PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat, MPI_Comm, MatReuse, Mat *);
805a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat, PetscInt, IS[], PetscInt);
815a576424SJed Brown PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat, PetscInt, IS *);
82b5b72c8aSIrina Sokolova PETSC_INTERN PetscErrorCode MatMPIBAIJSetPreallocation_MPIBAIJ(Mat B, PetscInt bs, PetscInt d_nz, const PetscInt *d_nnz, PetscInt o_nz, const PetscInt *o_nnz);
834de5dceeSHong Zhang PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIBAIJ(Mat, const PetscInt *, Mat, const PetscInt *, PetscInt *);
842726fb6dSPierre Jolivet 
852726fb6dSPierre Jolivet PETSC_INTERN PetscErrorCode MatConjugate_SeqBAIJ(Mat);
86