xref: /petsc/src/mat/impls/aij/seq/aij.h (revision 28b2fa4ae375c43e60412405d68670f082d116bf)
1b8a66259SBarry Smith 
22d40f771SBarry Smith #if !defined(__AIJ_H)
32d40f771SBarry Smith #define __AIJ_H
4e6b907acSBarry Smith 
57c4f633dSBarry Smith #include "private/matimpl.h"
6e6b907acSBarry Smith /*
7e6b907acSBarry Smith     Struct header shared by SeqAIJ, SeqBAIJ and SeqSBAIJ matrix formats
8e6b907acSBarry Smith */
9421e10b8SBarry Smith #define SEQAIJHEADER(datatype)	\
10e6b907acSBarry Smith   PetscTruth        roworiented;      /* if true, row-oriented input, default */\
11e6b907acSBarry Smith   PetscInt          nonew;            /* 1 don't add new nonzeros, -1 generate error on new */\
12*28b2fa4aSMatthew Knepley   PetscInt          nounused;         /* -1 generate error on unused space */\
13e6b907acSBarry Smith   PetscTruth        singlemalloc;     /* if true a, i, and j have been obtained with one big malloc */\
14e6b907acSBarry Smith   PetscInt          maxnz;            /* allocated nonzeros */\
15e6b907acSBarry Smith   PetscInt          *imax;            /* maximum space allocated for each row */\
16e6b907acSBarry Smith   PetscInt          *ilen;            /* actual length of each row */\
17e6b907acSBarry Smith   PetscInt          reallocs;         /* number of mallocs done during MatSetValues() \
18e6b907acSBarry Smith                                         as more values are set than were prealloced */\
19e6b907acSBarry Smith   PetscInt          rmax;             /* max nonzeros in any row */\
20e6b907acSBarry Smith   PetscTruth        keepzeroedrows;   /* keeps matrix structure same in calls to MatZeroRows()*/\
21e6b907acSBarry Smith   PetscTruth        ignorezeroentries; \
22e6b907acSBarry Smith   PetscInt          *xtoy,*xtoyB;     /* map nonzero pattern of X into Y's, used by MatAXPY() */\
23e6b907acSBarry Smith   Mat               XtoY;             /* used by MatAXPY() */\
24e6b907acSBarry Smith   PetscTruth        free_ij;          /* free the column indices j and row offsets i when the matrix is destroyed */ \
25e6b907acSBarry Smith   PetscTruth        free_a;           /* free the numerical values when matrix is destroy */ \
26e6b907acSBarry Smith   Mat_CompressedRow compressedrow;    /* use compressed row format */                      \
27e6b907acSBarry Smith   PetscInt          nz;               /* nonzeros */                                       \
28e6b907acSBarry Smith   PetscInt          *i;               /* pointer to beginning of each row */               \
29e6b907acSBarry Smith   PetscInt          *j;               /* column values: j + i[k] - 1 is start of row k */  \
30e6b907acSBarry Smith   PetscInt          *diag;            /* pointers to diagonal elements */                  \
31421e10b8SBarry Smith   datatype          *a;               /* nonzero elements */                               \
32e6b907acSBarry Smith   PetscScalar       *solve_work;      /* work space used in MatSolve */                    \
33e6b907acSBarry Smith   IS                row, col, icol    /* index sets, used for reorderings */
34e6b907acSBarry Smith 
352d40f771SBarry Smith /*
36ec8511deSBarry Smith   MATSEQAIJ format - Compressed row storage (also called Yale sparse matrix
37e6b907acSBarry Smith   format) or compressed sparse row (CSR).  The i[] and j[] arrays start at 0. For example,
38dfbc5765Svictorle   j[i[k]+p] is the pth column in row k.  Note that the diagonal
395768c4f9SLois Curfman McInnes   matrix elements are stored with the rest of the nonzeros (not separately).
402d40f771SBarry Smith */
41d35516d3SLois Curfman McInnes 
42e6b907acSBarry Smith /* Info about i-nodes (identical nodes) helper class for SeqAIJ */
43b8a66259SBarry Smith typedef struct {
44a77337e4SBarry Smith   MatScalar   *bdiag,*ibdiag;               /* diagonal blocks of matrix used for MatRelax_Inode() */
45f0d39aaaSBarry Smith   PetscInt    bdiagsize;                       /* length of bdiag and ibdiag */
4671f1c65dSBarry Smith   PetscTruth  ibdiagvalid;                     /* do ibdiag[] and bdiag[] contain the most recent values */
47f0d39aaaSBarry Smith 
48e6b907acSBarry Smith   PetscTruth use;
49e6b907acSBarry Smith   PetscInt   node_count;                    /* number of inodes */
50e6b907acSBarry Smith   PetscInt   *size;                         /* size of each inode */
51e6b907acSBarry Smith   PetscInt   limit;                         /* inode limit */
52e6b907acSBarry Smith   PetscInt   max_limit;                     /* maximum supported inode limit */
53e6b907acSBarry Smith   PetscTruth checked;                       /* if inodes have been checked for */
54e6b907acSBarry Smith } Mat_Inode;
55e6b907acSBarry Smith 
56e6b907acSBarry Smith EXTERN PetscErrorCode MatView_Inode(Mat,PetscViewer);
57e6b907acSBarry Smith EXTERN PetscErrorCode MatAssemblyEnd_Inode(Mat,MatAssemblyType);
58e6b907acSBarry Smith EXTERN PetscErrorCode MatDestroy_Inode(Mat);
59e6b907acSBarry Smith EXTERN PetscErrorCode MatCreate_Inode(Mat);
604e0d8c25SBarry Smith EXTERN PetscErrorCode MatSetOption_Inode(Mat,MatOption,PetscTruth);
61719d5645SBarry Smith EXTERN PetscErrorCode MatDuplicate_Inode(Mat,MatDuplicateOption,Mat*);
620481f469SBarry Smith EXTERN PetscErrorCode MatILUDTFactor_Inode(Mat,IS,IS,const MatFactorInfo*,Mat*);
630481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*);
640481f469SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_Inode(Mat,Mat,IS,IS,const MatFactorInfo*);
65e6b907acSBarry Smith 
66e6b907acSBarry Smith 
67e6b907acSBarry Smith typedef struct {
6854f21887SBarry Smith   SEQAIJHEADER(MatScalar);
69e6b907acSBarry Smith   Mat_Inode    inode;
7054f21887SBarry Smith   MatScalar    *saved_values;             /* location for stashing nonzero values of matrix */
7171f1c65dSBarry Smith 
7271f1c65dSBarry Smith   PetscScalar  *idiag,*mdiag,*ssor_work;  /* inverse of diagonal entries, diagonal values and workspace for Eisenstat trick */
7371f1c65dSBarry Smith   PetscTruth   idiagvalid;                     /* current idiag[] and mdiag[] are valid */
7471f1c65dSBarry Smith   PetscScalar  fshift,omega;                   /* last used omega and fshift */
7571f1c65dSBarry Smith 
763a7fca6bSBarry Smith   ISColoring   coloring;                  /* set with MatADSetColoring() used by MatADSetValues() */
77ec8511deSBarry Smith } Mat_SeqAIJ;
78b8a66259SBarry Smith 
79e6b907acSBarry Smith /*
80e6b907acSBarry Smith     Frees the a, i, and j arrays from the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
81e6b907acSBarry Smith */
82e6b907acSBarry Smith #undef __FUNCT__
83e6b907acSBarry Smith #define __FUNCT__ "MatSeqXAIJFreeAIJ"
84d0f46423SBarry Smith PETSC_STATIC_INLINE PetscErrorCode MatSeqXAIJFreeAIJ(Mat AA,MatScalar **a,PetscInt **j,PetscInt **i)
85d0f46423SBarry Smith {
86e6b907acSBarry Smith                                      PetscErrorCode ierr;
87e6b907acSBarry Smith                                      Mat_SeqAIJ     *A = (Mat_SeqAIJ*) AA->data;
88e6b907acSBarry Smith                                      if (A->singlemalloc) {
89e6b907acSBarry Smith                                        ierr = PetscFree3(*a,*j,*i);CHKERRQ(ierr);
90e6b907acSBarry Smith                                      } else {
91e6b907acSBarry Smith                                        if (A->free_a  && *a) {ierr = PetscFree(*a);CHKERRQ(ierr);}
92e6b907acSBarry Smith                                        if (A->free_ij && *j) {ierr = PetscFree(*j);CHKERRQ(ierr);}
93e6b907acSBarry Smith                                        if (A->free_ij && *i) {ierr = PetscFree(*i);CHKERRQ(ierr);}
94e6b907acSBarry Smith                                      }
95e6b907acSBarry Smith                                      *a = 0; *j = 0; *i = 0;
96e6b907acSBarry Smith                                      return 0;
97e6b907acSBarry Smith                                    }
98e6b907acSBarry Smith 
99e6b907acSBarry Smith /*
100e6b907acSBarry Smith     Allocates larger a, i, and j arrays for the XAIJ (AIJ, BAIJ, and SBAIJ) matrix types
101e6b907acSBarry Smith */
102421e10b8SBarry Smith #define MatSeqXAIJReallocateAIJ(Amat,AM,BS2,NROW,ROW,COL,RMAX,AA,AI,AJ,RP,AP,AIMAX,NONEW,datatype) \
103e6b907acSBarry Smith   if (NROW >= RMAX) {\
104e6b907acSBarry Smith 	Mat_SeqAIJ *Ain = (Mat_SeqAIJ*)Amat->data;\
105e6b907acSBarry Smith         /* there is no extra room in row, therefore enlarge */ \
106e6b907acSBarry Smith         PetscInt   new_nz = AI[AM] + CHUNKSIZE,len,*new_i=0,*new_j=0; \
107421e10b8SBarry Smith         datatype   *new_a; \
108e6b907acSBarry Smith  \
109e6b907acSBarry Smith         if (NONEW == -2) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"New nonzero at (%D,%D) caused a malloc",ROW,COL); \
110e6b907acSBarry Smith         /* malloc new storage space */ \
111421e10b8SBarry Smith         ierr = PetscMalloc3(BS2*new_nz,datatype,&new_a,new_nz,PetscInt,&new_j,AM+1,PetscInt,&new_i);CHKERRQ(ierr);\
112e6b907acSBarry Smith  \
113e6b907acSBarry Smith         /* copy over old data into new slots */ \
114e6b907acSBarry Smith         for (ii=0; ii<ROW+1; ii++) {new_i[ii] = AI[ii];} \
115e6b907acSBarry Smith         for (ii=ROW+1; ii<AM+1; ii++) {new_i[ii] = AI[ii]+CHUNKSIZE;} \
116e6b907acSBarry Smith         ierr = PetscMemcpy(new_j,AJ,(AI[ROW]+NROW)*sizeof(PetscInt));CHKERRQ(ierr); \
117e6b907acSBarry Smith         len = (new_nz - CHUNKSIZE - AI[ROW] - NROW); \
118e6b907acSBarry Smith         ierr = PetscMemcpy(new_j+AI[ROW]+NROW+CHUNKSIZE,AJ+AI[ROW]+NROW,len*sizeof(PetscInt));CHKERRQ(ierr); \
119421e10b8SBarry Smith         ierr = PetscMemcpy(new_a,AA,BS2*(AI[ROW]+NROW)*sizeof(datatype));CHKERRQ(ierr); \
120421e10b8SBarry Smith         ierr = PetscMemzero(new_a+BS2*(AI[ROW]+NROW),BS2*CHUNKSIZE*sizeof(datatype));CHKERRQ(ierr);\
121421e10b8SBarry Smith         ierr = PetscMemcpy(new_a+BS2*(AI[ROW]+NROW+CHUNKSIZE),AA+BS2*(AI[ROW]+NROW),BS2*len*sizeof(datatype));CHKERRQ(ierr);  \
122e6b907acSBarry Smith         /* free up old matrix storage */ \
123e6b907acSBarry Smith         ierr = MatSeqXAIJFreeAIJ(A,&Ain->a,&Ain->j,&Ain->i);CHKERRQ(ierr);\
124421e10b8SBarry Smith         AA = new_a; \
125421e10b8SBarry Smith         Ain->a = (MatScalar*) new_a;		   \
126421e10b8SBarry Smith         AI = Ain->i = new_i; AJ = Ain->j = new_j;  \
127e6b907acSBarry Smith         Ain->singlemalloc = PETSC_TRUE; \
128e6b907acSBarry Smith  \
129e6b907acSBarry Smith         RP          = AJ + AI[ROW]; AP = AA + BS2*AI[ROW]; \
130e6b907acSBarry Smith         RMAX        = AIMAX[ROW] = AIMAX[ROW] + CHUNKSIZE; \
131e6b907acSBarry Smith         Ain->maxnz += CHUNKSIZE; \
132e6b907acSBarry Smith         Ain->reallocs++; \
133e6b907acSBarry Smith       } \
134e6b907acSBarry Smith 
135f349cc81SSatish Balay EXTERN_C_BEGIN
136ab93d7beSBarry Smith EXTERN PetscErrorCode MatSeqAIJSetPreallocation_SeqAIJ(Mat,PetscInt,PetscInt*);
137f349cc81SSatish Balay EXTERN_C_END
1380481f469SBarry Smith EXTERN PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
1390481f469SBarry Smith EXTERN PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
1400481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat,Mat,IS,const MatFactorInfo*);
1410481f469SBarry Smith EXTERN PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
142dfbe8321SBarry Smith EXTERN PetscErrorCode MatDuplicate_SeqAIJ(Mat,MatDuplicateOption,Mat*);
14309f38230SBarry Smith EXTERN PetscErrorCode MatMissingDiagonal_SeqAIJ(Mat,PetscTruth*,PetscInt*);
144dfbe8321SBarry Smith EXTERN PetscErrorCode MatMarkDiagonal_SeqAIJ(Mat);
14508480c60SBarry Smith 
146dfbe8321SBarry Smith EXTERN PetscErrorCode MatMult_SeqAIJ(Mat A,Vec,Vec);
147dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
148dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTranspose_SeqAIJ(Mat A,Vec,Vec);
149dfbe8321SBarry Smith EXTERN PetscErrorCode MatMultTransposeAdd_SeqAIJ(Mat A,Vec,Vec,Vec);
15038baddfdSBarry Smith EXTERN PetscErrorCode MatRelax_SeqAIJ(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
15108480c60SBarry Smith 
152dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetColoring_SeqAIJ(Mat,ISColoring);
153dfbe8321SBarry Smith EXTERN PetscErrorCode MatSetValuesAdic_SeqAIJ(Mat,void*);
15438baddfdSBarry Smith EXTERN PetscErrorCode MatSetValuesAdifor_SeqAIJ(Mat,PetscInt,void*);
1553a7fca6bSBarry Smith 
15638baddfdSBarry Smith EXTERN PetscErrorCode MatGetSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
1570390132cSHong Zhang EXTERN PetscErrorCode MatGetSymbolicTransposeReduced_SeqAIJ(Mat,PetscInt,PetscInt,PetscInt *[],PetscInt *[]);
15838baddfdSBarry Smith EXTERN PetscErrorCode MatRestoreSymbolicTranspose_SeqAIJ(Mat,PetscInt *[],PetscInt *[]);
15938baddfdSBarry Smith EXTERN PetscErrorCode MatToSymmetricIJ_SeqAIJ(PetscInt,PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt**,PetscInt**);
1600481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat,Mat,IS,IS,const MatFactorInfo*);
1610481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat,Mat,const MatFactorInfo*);
1620481f469SBarry Smith EXTERN PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat,Mat,const MatFactorInfo*);
1630481f469SBarry Smith EXTERN PetscErrorCode MatLUFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*);
164dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ(Mat,Vec,Vec);
1658d8c7c43SBarry Smith EXTERN PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat,Vec,Vec);
166137fb511SHong Zhang EXTERN PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat,Vec,Vec);
167dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveAdd_SeqAIJ(Mat,Vec,Vec,Vec);
168dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTranspose_SeqAIJ(Mat,Vec,Vec);
169dfbe8321SBarry Smith EXTERN PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat,Vec,Vec,Vec);
1702bebee5dSHong Zhang EXTERN PetscErrorCode MatMatSolve_SeqAIJ(Mat,Mat,Mat);
171dfbe8321SBarry Smith EXTERN PetscErrorCode MatEqual_SeqAIJ(Mat A,Mat B,PetscTruth* flg);
172dfbe8321SBarry Smith EXTERN PetscErrorCode MatFDColoringCreate_SeqAIJ(Mat,ISColoring,MatFDColoring);
1730481f469SBarry Smith EXTERN PetscErrorCode MatILUDTFactor_SeqAIJ(Mat,IS,IS,const MatFactorInfo*,Mat*);
174a313700dSBarry Smith EXTERN PetscErrorCode MatLoad_SeqAIJ(PetscViewer, const MatType,Mat*);
175dfbe8321SBarry Smith EXTERN PetscErrorCode RegisterApplyPtAPRoutines_Private(Mat);
176dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
177dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
178dfbe8321SBarry Smith EXTERN PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
1797ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ(Mat,Mat,PetscReal,Mat*);
1807ba1a0bfSKris Buschelman EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ(Mat,Mat,Mat);
181dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
182dfbe8321SBarry Smith EXTERN PetscErrorCode MatPtAPNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
183bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTranspose_SeqAIJ_SeqAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
184bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeSymbolic_SeqAIJ_SeqAIJ(Mat,Mat,PetscReal,Mat*);
185bc011b1eSHong Zhang EXTERN PetscErrorCode MatMatMultTransposeNumeric_SeqAIJ_SeqAIJ(Mat,Mat,Mat);
18638baddfdSBarry Smith EXTERN PetscErrorCode MatSetValues_SeqAIJ(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
187a77337e4SBarry Smith EXTERN PetscErrorCode MatGetRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
188a77337e4SBarry Smith EXTERN PetscErrorCode MatRestoreRow_SeqAIJ(Mat,PetscInt,PetscInt*,PetscInt**,PetscScalar**);
189f4df32b1SMatthew Knepley EXTERN PetscErrorCode MatAXPY_SeqAIJ(Mat,PetscScalar,Mat,MatStructure);
1908f7157efSSatish Balay EXTERN PetscErrorCode MatGetRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
1918f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreRowIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
1928f7157efSSatish Balay EXTERN PetscErrorCode MatGetColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
1938f7157efSSatish Balay EXTERN PetscErrorCode MatRestoreColumnIJ_SeqAIJ(Mat,PetscInt,PetscTruth,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
194b24902e0SBarry Smith EXTERN PetscErrorCode MatDestroy_SeqAIJ(Mat);
195b24902e0SBarry Smith EXTERN PetscErrorCode MatView_SeqAIJ(Mat,PetscViewer);
1969af31e4aSHong Zhang 
19797304618SKris Buschelman EXTERN_C_BEGIN
1988cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqSBAIJ(Mat,const MatType,MatReuse,Mat*);
1998cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqBAIJ(Mat,const MatType,MatReuse,Mat*);
2008cf70c4bSSatish Balay EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatConvert_SeqAIJ_SeqCSRPERM(Mat,const MatType,MatReuse,Mat*);
201be1d678aSKris Buschelman EXTERN PetscErrorCode PETSCMAT_DLLEXPORT MatReorderForNonzeroDiagonal_SeqAIJ(Mat,PetscReal,IS,IS);
20297304618SKris Buschelman EXTERN_C_END
20370f19b1fSKris Buschelman 
2042d40f771SBarry Smith #endif
205