xref: /petsc/src/mat/impls/aij/seq/aij.h (revision d32f9abdbc052d6e1fd06679b17a55415c3aae30)
1 
2 #if !defined(__AIJ_H)
3 #define __AIJ_H
4 
5 #include "include/private/matimpl.h"
6 /*
7     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
8 */
9 #define SEQAIJHEADER(datatype)	\
10   PetscTruth        roworiented;      /* if true, row-oriented input, default */\
11   PetscInt          nonew;            /* 1 don't add new nonzeros, -1 generate error on new */\
12   PetscTruth        singlemalloc;     /* if true a, i, and j have been obtained with one big malloc */\
13   PetscInt          maxnz;            /* allocated nonzeros */\
14   PetscInt          *imax;            /* maximum space allocated for each row */\
15   PetscInt          *ilen;            /* actual length of each row */\
16   PetscInt          reallocs;         /* number of mallocs done during MatSetValues() \
17                                         as more values are set than were prealloced */\
18   PetscInt          rmax;             /* max nonzeros in any row */\
19   PetscTruth        keepzeroedrows;   /* keeps matrix structure same in calls to MatZeroRows()*/\
20   PetscTruth        ignorezeroentries; \
21   PetscInt          *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */\
22   Mat               XtoY;             /* used by MatAXPY() */\
23   PetscTruth        free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
24   PetscTruth        free_a;           /* free the numerical values when matrix is destroy */ \
25   Mat_CompressedRow compressedrow;    /* use compressed row format */                      \
26   PetscInt          nz;               /* nonzeros */                                       \
27   PetscInt          *i;               /* pointer to beginning of each row */               \
28   PetscInt          *j;               /* column values: j + i[k] - 1 is start of row k */  \
29   PetscInt          *diag;            /* pointers to diagonal elements */                  \
30   datatype          *a;               /* nonzero elements */                               \
31   PetscScalar       *solve_work;      /* work space used in MatSolve */                    \
32   IS                row, col, icol    /* index sets, used for reorderings */
33 
34 /*
35   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
36   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
37   j[i[k]+p] is the pth column in row k.  Note that the diagonal
38   matrix elements are stored with the rest of the nonzeros (not separately).
39 */
40 
41 /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
42 typedef struct {
43   MatScalar   *bdiag,*ibdiag;               /* diagonal blocks of matrix used for MatRelax_Inode() */
44   PetscInt    bdiagsize;                       /* length of bdiag and ibdiag */
45   PetscTruth  ibdiagvalid;                     /* do ibdiag[] and bdiag[] contain the most recent values */
46 
47   PetscTruth use;
48   PetscInt   node_count;                    /* number of inodes */
49   PetscInt   *size;                         /* size of each inode */
50   PetscInt   limit;                         /* inode limit */
51   PetscInt   max_limit;                     /* maximum supported inode limit */
52   PetscTruth checked;                       /* if inodes have been checked for */
53 } Mat_Inode;
54 
55 EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer);
56 EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType);
57 EXTERN PetscErrorCode MatDestroy_Inode(Mat);
58 EXTERN PetscErrorCode MatCreate_Inode(Mat);
59 EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth);
60 EXTERN PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *B);
61 EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
62 EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
63 EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
64 
65 
66 typedef struct {
67   SEQAIJHEADER(MatScalar);
68   Mat_Inode    inode;
69   MatScalar    *saved_values;             /* location for stashing nonzero values of matrix */
70 
71   PetscScalar  *idiag,*mdiag,*ssor_work;  /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
72   PetscTruth   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
73   PetscScalar  fshift,omega;                   /* last used omega and fshift */
74 
75   ISColoring   coloring;                  /* set with MatADSetColoring() used by MatADSetValues() */
76 } Mat_SeqAIJ;
77 
78 /*
79     Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
80 */
81 #undef __FUNCT__
82 #define __FUNCT__ "MatSeqXAIJFreeAIJ"
83 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i) {
84                                      PetscErrorCode ierr;
85                                      Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
86                                      if (A->singlemalloc) {
87                                        ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
88                                      } else {
89                                        if (A->free_a  && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);}
90                                        if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);}
91                                        if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);}
92                                      }
93                                      *a = 0; *j = 0; *i = 0;
94                                      return 0;
95                                    }
96 
97 /*
98     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
99 */
100 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
101   if (NROW >= RMAX) {\
102 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
103         /* there is no extra room in row, therefore enlarge */ \
104         PetscInt   new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
105         datatype   *new_a; \
106  \
107         if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
108         /* malloc new storage space */ \
109         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
110  \
111         /* copy over old data into new slots */ \
112         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
113         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
114         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
115         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
116         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
117         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
118         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
119         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
120         /* free up old matrix storage */ \
121         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
122         AA = new_a; \
123         Ain->a = (MatScalar*) new_a;		   \
124         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
125         Ain->singlemalloc = PETSC_TRUE; \
126  \
127         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
128         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
129         Ain->maxnz += CHUNKSIZE; \
130         Ain->reallocs++; \
131       } \
132 
133 EXTERN_C_BEGIN
134 EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*);
135 EXTERN_C_END
136 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
137 EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat *);
138 EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat*);
139 EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*);
140 EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
141 EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*);
142 EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
143 
144 EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
145 EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
146 EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
147 EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
148 EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
149 
150 EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
151 EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
152 EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
153 
154 EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
155 EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
156 EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
157 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
158 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
159 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*);
160 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,MatFactorInfo*,Mat*);
161 EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*);
162 EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
163 EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
164 EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
165 EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
166 EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
167 EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
168 EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
169 EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
170 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
171 EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*);
172 EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
173 EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
174 EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
175 EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
176 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
177 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
178 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
179 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
180 EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
181 EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
182 EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
183 EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
184 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
185 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
186 EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
187 EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
188 EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
189 EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
190 EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
191 EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
192 EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
193 
194 EXTERN_C_BEGIN
195 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*);
196 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*);
197 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*);
198 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
199 EXTERN_C_END
200 
201 #endif
202