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