Actual source code: mpibaij.h
4: #include src/mat/impls/baij/seq/baij.h
5: #include src/sys/ctable.h
7: #if defined (PETSC_USE_CTABLE)
8: #define PETSCTABLE PetscTable
9: #else
10: #define PETSCTABLE PetscInt*
11: #endif
13: #define MPIBAIJHEADER \
14: PetscInt *rangebs; /* rmap.range/bs */ \
15: PetscInt rstartbs,rendbs,cstartbs,cendbs; /* map values / bs */ \
16: Mat A,B; /* local submatrices: A (diag part), B (off-diag part) */ \
17: PetscMPIInt size; /* size of communicator */ \
18: PetscMPIInt rank; /* rank of proc in communicator */ \
19: PetscInt bs2; /* block size, bs2 = bs*bs */ \
20: PetscInt Mbs,Nbs; /* number block rows/cols in matrix; M/bs, N/bs */ \
21: PetscInt mbs,nbs; /* number block rows/cols on processor; m/bs, n/bs */ \
22: \
23: /* The following variables are used for matrix assembly */ \
24: \
25: PetscTruth donotstash; /* if 1, off processor entries dropped */ \
26: MPI_Request *send_waits; /* array of send requests */ \
27: MPI_Request *recv_waits; /* array of receive requests */ \
28: PetscInt nsends,nrecvs; /* numbers of sends and receives */ \
29: MatScalar *svalues,*rvalues; /* sending and receiving data */ \
30: PetscInt rmax; /* maximum message length */ \
31: PETSCTABLE colmap; /* local col number of off-diag col */ \
32: \
33: PetscInt *garray; /* work array */ \
34: \
35: /* The following variable is used by blocked matrix assembly */ \
36: MatScalar *barray; /* Block array of size bs2 */ \
37: \
38: /* The following variables are used for matrix-vector products */ \
39: \
40: Vec lvec; /* local vector */ \
41: VecScatter Mvctx; /* scatter context for vector */ \
42: PetscTruth roworiented; /* if true, row-oriented input, default true */ \
43: \
44: /* The following variables are for MatGetRow() */ \
45: \
46: PetscInt *rowindices; /* column indices for row */ \
47: PetscScalar *rowvalues; /* nonzero values in row */ \
48: PetscTruth getrowactive; /* indicates MatGetRow(), not restored */ \
49: \
50: /* Some variables to make MatSetValues and others more efficient */ \
51: PetscInt rstart_bs,rend_bs; \
52: PetscInt cstart_bs,cend_bs; \
53: PetscInt *ht; /* Hash table to speed up matrix assembly */ \
54: MatScalar **hd; /* Hash table data */ \
55: PetscInt ht_size; \
56: PetscInt ht_total_ct,ht_insert_ct; /* Hash table statistics */ \
57: PetscTruth ht_flag; /* Flag to indicate if hash tables are used */ \
58: double ht_fact; /* Factor to determine the HT size */ \
59: \
60: PetscInt setvalueslen; /* only used for single precision computations */ \
61: MatScalar *setvaluescopy /* area double precision values in MatSetValuesXXX() are copied*/ \
62: /* before calling MatSetValuesXXX_MPIBAIJ_MatScalar() */
64: typedef struct {
65: MPIBAIJHEADER;
66: } Mat_MPIBAIJ;
68: EXTERN PetscErrorCode MatLoad_MPIBAIJ(PetscViewer, MatType,Mat*);
69: EXTERN PetscErrorCode CreateColmap_MPIBAIJ_Private(Mat);
70: EXTERN PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
71: #endif