Actual source code: mpibaij.h

  1: #pragma once
  2: #include <../src/mat/impls/baij/seq/baij.h>
  3: #include <../src/mat/impls/aij/mpi/mpiaij.h>
  4: #include <petsc/private/hashmapi.h>

  6: #define MPIBAIJHEADER \
  7:   MPIAIJHEADER; \
  8:   PetscInt *rangebs;                            /* rmap->range/bs */ \
  9:   PetscInt  rstartbs, rendbs, cstartbs, cendbs; /* map values / bs  */ \
 10:   PetscInt  bs2;                                /* block size, bs2 = bs*bs */ \
 11:   PetscInt  Mbs, Nbs;                           /* number block rows/cols in matrix; M/bs, N/bs */ \
 12:   PetscInt  mbs, nbs;                           /* number block rows/cols on processor; m/bs, n/bs */ \
 13: \
 14:   /* The following variables are used for matrix assembly */ \
 15:   PetscBool subset_off_proc_entries; /* PETSC_TRUE if assembly will always communicate a subset of the entries communicated the first time */ \
 16: \
 17:   /* The following variable is used by blocked matrix assembly */ \
 18:   MatScalar *barray; /* Block array of size bs2 */ \
 19: \
 20:   /* Some variables to make MatSetValues and others more efficient */ \
 21:   PetscInt    rstart_bs, rend_bs; \
 22:   PetscInt    cstart_bs, cend_bs; \
 23:   PetscInt   *ht; /* Hash table to speed up matrix assembly */ \
 24:   MatScalar **hd; /* Hash table data */ \
 25:   PetscInt    ht_size; \
 26:   PetscInt    ht_total_ct, ht_insert_ct; /* Hash table statistics */ \
 27:   PetscBool   ht_flag;                   /* Flag to indicate if hash tables are used */ \
 28:   double      ht_fact;                   /* Factor to determine the HT size */ \
 29: \
 30:   PetscInt   setvalueslen;  /* only used for single precision computations */ \
 31:   MatScalar *setvaluescopy; /* area double precision values in MatSetValuesXXX() are copied*/ \
 32:                             /* before calling MatSetValuesXXX_MPIBAIJ_MatScalar() */ \
 33:   PetscBool      ijonly;    /* used in  MatCreateSubMatrices_MPIBAIJ_local() for getting ij structure only */ \
 34:   struct _MatOps cops

 36: typedef struct {
 37:   MPIBAIJHEADER;
 38: } Mat_MPIBAIJ;

 40: PETSC_INTERN PetscErrorCode MatView_MPIBAIJ(Mat, PetscViewer);
 41: PETSC_INTERN PetscErrorCode MatLoad_MPIBAIJ(Mat, PetscViewer);
 42: PETSC_INTERN PetscErrorCode MatView_MPIBAIJ_Binary(Mat, PetscViewer);
 43: PETSC_INTERN PetscErrorCode MatLoad_MPIBAIJ_Binary(Mat, PetscViewer);

 45: PETSC_INTERN PetscErrorCode MatCreateColmap_MPIBAIJ_Private(Mat);
 46: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIBAIJ(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
 47: PETSC_INTERN PetscErrorCode MatCreateSubMatrices_MPIBAIJ_local(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *, PetscBool);
 48: PETSC_INTERN PetscErrorCode MatCreateSubMatrix_MPIBAIJ_Private(Mat, IS, IS, PetscInt, MatReuse, Mat *, PetscBool);
 49: PETSC_INTERN PetscErrorCode MatGetMultiProcBlock_MPIBAIJ(Mat, MPI_Comm, MatReuse, Mat *);
 50: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIBAIJ(Mat, PetscInt, IS[], PetscInt);
 51: PETSC_INTERN PetscErrorCode MatIncreaseOverlap_MPIBAIJ_Once(Mat, PetscInt, IS *);
 52: PETSC_INTERN PetscErrorCode MatMPIBAIJSetPreallocation_MPIBAIJ(Mat B, PetscInt bs, PetscInt d_nz, const PetscInt *d_nnz, PetscInt o_nz, const PetscInt *o_nnz);
 53: PETSC_INTERN PetscErrorCode MatAXPYGetPreallocation_MPIBAIJ(Mat, const PetscInt *, Mat, const PetscInt *, PetscInt *);

 55: PETSC_INTERN PetscErrorCode MatConjugate_SeqBAIJ(Mat);