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 *rowners,*cowners; /* ranges owned by each processor, in blocks */ \
15: PetscInt *rowners_bs; /* rowners*bs */ \
16: PetscInt rstart,rend; /* starting and ending owned rows */ \
17: PetscInt cstart,cend; /* starting and ending owned columns */ \
18: Mat A,B; /* local submatrices: A (diag part), B (off-diag part) */ \
19: PetscMPIInt size; /* size of communicator */ \
20: PetscMPIInt rank; /* rank of proc in communicator */ \
21: PetscInt bs2; /* block size, bs2 = bs*bs */ \
22: PetscInt Mbs,Nbs; /* number block rows/cols in matrix; M/bs, N/bs */ \
23: PetscInt mbs,nbs; /* number block rows/cols on processor; m/bs, n/bs */ \
24: \
25: /* The following variables are used for matrix assembly */ \
26: \
27: PetscTruth donotstash; /* if 1, off processor entries dropped */ \
28: MPI_Request *send_waits; /* array of send requests */ \
29: MPI_Request *recv_waits; /* array of receive requests */ \
30: PetscInt nsends,nrecvs; /* numbers of sends and receives */ \
31: MatScalar *svalues,*rvalues; /* sending and receiving data */ \
32: PetscInt rmax; /* maximum message length */ \
33: PETSCTABLE colmap; /* local col number of off-diag col */ \
34: \
35: PetscInt *garray; /* work array */ \
36: \
37: /* The following variable is used by blocked matrix assembly */ \
38: MatScalar *barray; /* Block array of size bs2 */ \
39: \
40: /* The following variables are used for matrix-vector products */ \
41: \
42: Vec lvec; /* local vector */ \
43: VecScatter Mvctx; /* scatter context for vector */ \
44: PetscTruth roworiented; /* if true, row-oriented input, default true */ \
45: \
46: /* The following variables are for MatGetRow() */ \
47: \
48: PetscInt *rowindices; /* column indices for row */ \
49: PetscScalar *rowvalues; /* nonzero values in row */ \
50: PetscTruth getrowactive; /* indicates MatGetRow(), not restored */ \
51: \
52: /* Some variables to make MatSetValues and others more efficient */ \
53: PetscInt rstart_bs,rend_bs; \
54: PetscInt cstart_bs,cend_bs; \
55: PetscInt *ht; /* Hash table to speed up matrix assembly */ \
56: MatScalar **hd; /* Hash table data */ \
57: PetscInt ht_size; \
58: PetscInt ht_total_ct,ht_insert_ct; /* Hash table statistics */ \
59: PetscTruth ht_flag; /* Flag to indicate if hash tables are used */ \
60: double ht_fact; /* Factor to determine the HT size */ \
61: \
62: PetscInt setvalueslen; /* only used for single precision computations */ \
63: MatScalar *setvaluescopy /* area double precision values in MatSetValuesXXX() are copied*/ \
64: /* before calling MatSetValuesXXX_MPIBAIJ_MatScalar() */
66: typedef struct {
67: MPIBAIJHEADER;
68: } Mat_MPIBAIJ;
70: EXTERN PetscErrorCode MatLoad_MPIBAIJ(PetscViewer, MatType,Mat*);
71: EXTERN PetscErrorCode CreateColmap_MPIBAIJ_Private(Mat);
72: EXTERN PetscErrorCode MatGetSubMatrices_MPIBAIJ(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat*[]);
73: #endif