xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 22612f2f7cceb60caedd65384cdf99fc989f2aeb)
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   PetscTruth use;
44   PetscInt   node_count;                    /* number of inodes */
45   PetscInt   *size;                         /* size of each inode */
46   PetscInt   limit;                         /* inode limit */
47   PetscInt   max_limit;                     /* maximum supported inode limit */
48   PetscTruth checked;                       /* if inodes have been checked for */
49 } Mat_Inode;
50 
51 EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer);
52 EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType);
53 EXTERN PetscErrorCode MatDestroy_Inode(Mat);
54 EXTERN PetscErrorCode MatCreate_Inode(Mat);
55 EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth);
56 EXTERN PetscErrorCode MatDuplicate_Inode(Mat A,MatDuplicateOption cpvalues,Mat *B);
57 EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
58 EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
59 EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat A,IS isrow,IS iscol,MatFactorInfo *info,Mat *fact);
60 
61 
62 typedef struct {
63   SEQAIJHEADER(PetscScalar);
64   Mat_Inode    inode;
65   PetscScalar  *saved_values;    /* location for stashing nonzero values of matrix */
66   PetscScalar  *idiag,*ssor;     /* inverse of diagonal entries; space for eisen */
67   ISColoring   coloring;         /* set with MatADSetColoring() used by MatADSetValues() */
68 } Mat_SeqAIJ;
69 
70 /*
71     Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
72 */
73 #undef __FUNCT__
74 #define __FUNCT__ "MatSeqXAIJFreeAIJ"
75 PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,PetscScalar **a,PetscInt **j,PetscInt **i) {
76                                      PetscErrorCode ierr;
77                                      Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
78                                      if (A->singlemalloc) {
79                                        ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
80                                      } else {
81                                        if (A->free_a  && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);}
82                                        if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);}
83                                        if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);}
84                                      }
85                                      *a = 0; *j = 0; *i = 0;
86                                      return 0;
87                                    }
88 
89 /*
90     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
91 */
92 #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
93   if (NROW >= RMAX) {\
94 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
95         /* there is no extra room in row, therefore enlarge */ \
96         PetscInt   new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
97         datatype   *new_a; \
98  \
99         if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
100         /* malloc new storage space */ \
101         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
102  \
103         /* copy over old data into new slots */ \
104         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
105         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
106         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
107         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
108         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
109         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
110         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
111         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
112         /* free up old matrix storage */ \
113         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
114         AA = new_a; \
115         Ain->a = (MatScalar*) new_a;		   \
116         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
117         Ain->singlemalloc = PETSC_TRUE; \
118  \
119         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
120         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
121         Ain->maxnz += CHUNKSIZE; \
122         Ain->reallocs++; \
123       } \
124 
125 EXTERN_C_BEGIN
126 EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*);
127 EXTERN_C_END
128 EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat *);
129 EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat *);
130 EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,IS,MatFactorInfo*,Mat*);
131 EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*);
132 EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
133 EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*);
134 EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
135 
136 EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
137 EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
138 EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
139 EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
140 EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
141 
142 EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
143 EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
144 EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
145 
146 EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
147 EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
148 EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
149 EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
150 EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
151 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,MatFactorInfo*,Mat*);
152 EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,MatFactorInfo*,Mat*);
153 EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*);
154 EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
155 EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
156 EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
157 EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
158 EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
159 EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
160 EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
161 EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
162 EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,MatFactorInfo*,Mat*);
163 EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, MatType,Mat*);
164 EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
165 EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
166 EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
167 EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
168 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
169 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
170 EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
171 EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
172 EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
173 EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
174 EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
175 EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
176 EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
177 EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
178 EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
179 EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
180 EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
181 EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
182 EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
183 
184 EXTERN_C_BEGIN
185 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat, MatType,MatReuse,Mat*);
186 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat, MatType,MatReuse,Mat*);
187 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,MatType,MatReuse,Mat*);
188 EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
189 EXTERN_C_END
190 
191 #endif
192