xref: /petsc/src/mat/impls/aij/seq/bas/spbas.h (revision 9767453c85a88b60d52ea72032faaafc609f88f5)
1b2f1ab58SBarry Smith /*
2b2f1ab58SBarry Smith    Define type spbas_matrix: sparse matrices using pointers
3b2f1ab58SBarry Smith 
4b2f1ab58SBarry Smith    Global matrix information
5b2f1ab58SBarry Smith       nrows, ncols: dimensions
6b2f1ab58SBarry Smith       nnz         : number of nonzeros (in whole matrix)
7b2f1ab58SBarry Smith       col_idx_type: storage scheme for column numbers
8b2f1ab58SBarry Smith                     SPBAS_COLUMN_NUMBERS:
9b2f1ab58SBarry Smith                         array icol contains column indices:
10b2f1ab58SBarry Smith                            A(i,icol[i][j]) = values[i][j]
11b2f1ab58SBarry Smith                     SPBAS_DIAGONAL_OFFSETS:
12b2f1ab58SBarry Smith                         array icol contains diagonal offsets:
13b2f1ab58SBarry Smith                            A(i,i+icol[i][j]) = values[i][j]
14b2f1ab58SBarry Smith                     SPBAS_OFFSET_ARRAY:
15b2f1ab58SBarry Smith                         array icol contains offsets wrt array
16b2f1ab58SBarry Smith                         icol0:
17b2f1ab58SBarry Smith                            A(i,icol0[i]+icol[i][j]) = values[i][j]
18b2f1ab58SBarry Smith 
19b2f1ab58SBarry Smith    Information about each row
20b2f1ab58SBarry Smith       row_nnz     : number of nonzeros for each row
21b2f1ab58SBarry Smith       icol0       : column index offset (when needed, otherwise NULL)
22b2f1ab58SBarry Smith       icols       : array of diagonal offsets for each row, as descibed
23b2f1ab58SBarry Smith                     for col_idx_type, above
24b2f1ab58SBarry Smith       values      : array of matrix entries for each row
25b2f1ab58SBarry Smith                     when values == NULL, this matrix is really
26b2f1ab58SBarry Smith                     a sparseness pattern, not a matrix
27b2f1ab58SBarry Smith 
28b2f1ab58SBarry Smith    The other fields describe the way in which the data are stored
29b2f1ab58SBarry Smith    in memory.
30b2f1ab58SBarry Smith 
31b2f1ab58SBarry Smith       block_data  : The pointers icols[i] all point to places in a
32b2f1ab58SBarry Smith                     single allocated array. Only for icols[0] was
33b2f1ab58SBarry Smith                     malloc called. Freeing icols[0] will free
34b2f1ab58SBarry Smith                     all other icols=arrays as well.
35b2f1ab58SBarry Smith                     Same for arrays values[i]
36b2f1ab58SBarry Smith */
37b2f1ab58SBarry Smith 
38b2f1ab58SBarry Smith #define SPBAS_COLUMN_NUMBERS   (0)
39b2f1ab58SBarry Smith #define SPBAS_DIAGONAL_OFFSETS (1)
40b2f1ab58SBarry Smith #define SPBAS_OFFSET_ARRAY     (2)
41b2f1ab58SBarry Smith 
42b2f1ab58SBarry Smith #define NEGATIVE_DIAGONAL (-42)
43b2f1ab58SBarry Smith 
44b2f1ab58SBarry Smith typedef struct {
45b2f1ab58SBarry Smith  PetscInt    nrows;
46b2f1ab58SBarry Smith  PetscInt    ncols;
47b2f1ab58SBarry Smith  PetscInt    nnz;
48b2f1ab58SBarry Smith  PetscInt    col_idx_type;
49b2f1ab58SBarry Smith 
50b2f1ab58SBarry Smith  PetscInt    *row_nnz;
51b2f1ab58SBarry Smith  PetscInt    *icol0;
52b2f1ab58SBarry Smith  PetscInt    **icols;
53b2f1ab58SBarry Smith  PetscScalar **values;
54b2f1ab58SBarry Smith 
55b2f1ab58SBarry Smith  PetscTruth  block_data;
56b2f1ab58SBarry Smith  PetscInt    n_alloc_icol;
57b2f1ab58SBarry Smith  PetscInt    n_alloc_val;
58b2f1ab58SBarry Smith  PetscInt    *alloc_icol;
59b2f1ab58SBarry Smith  PetscScalar *alloc_val;
60b2f1ab58SBarry Smith } spbas_matrix ;
61b2f1ab58SBarry Smith 
62b2f1ab58SBarry Smith 
63b2f1ab58SBarry Smith /*
64b2f1ab58SBarry Smith   spbas_compress_pattern:
65b2f1ab58SBarry Smith      calculate a compressed sparseness pattern for a sparseness pattern
66b2f1ab58SBarry Smith      given in compressed row storage. The compressed sparseness pattern may
67b2f1ab58SBarry Smith      require (much) less memory.
68b2f1ab58SBarry Smith 
69b2f1ab58SBarry Smith   spbas_memory_requirement:
70b2f1ab58SBarry Smith      Calculate the number of bytes needed to store tha matrix
71b2f1ab58SBarry Smith 
72b2f1ab58SBarry Smith   spbas_incomplete_cholesky:
73b2f1ab58SBarry Smith      Incomplete Cholesky decomposition
74b2f1ab58SBarry Smith 
75b2f1ab58SBarry Smith   spbas_delete:
76b2f1ab58SBarry Smith      de-allocate the arrays owned by this matrix
77b2f1ab58SBarry Smith 
78b2f1ab58SBarry Smith   spbas_matrix_to_crs:
79b2f1ab58SBarry Smith      Convert an spbas_matrix to compessed row storage
80b2f1ab58SBarry Smith 
81b2f1ab58SBarry Smith   spbas_dump:
82b2f1ab58SBarry Smith      Print the matrix in i,j,val-format
83b2f1ab58SBarry Smith 
84b2f1ab58SBarry Smith   spbas_transpose:
85b2f1ab58SBarry Smith      Return the transpose of a matrix
86b2f1ab58SBarry Smith 
87b2f1ab58SBarry Smith   spbas_pattern_only:
88b2f1ab58SBarry Smith      Return the sparseness pattern (matrix without values) of a
89b2f1ab58SBarry Smith      compressed row storage
90b2f1ab58SBarry Smith */
91*9767453cSBarry Smith PetscErrorCode spbas_compress_pattern(PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,spbas_matrix*,PetscReal*);
92b2f1ab58SBarry Smith long int spbas_memory_requirement( spbas_matrix);
93b2f1ab58SBarry Smith PetscErrorCode spbas_delete(spbas_matrix);
94d36aa34eSBarry Smith PetscErrorCode spbas_incomplete_cholesky(Mat,const PetscInt*,const PetscInt*,spbas_matrix,PetscReal,PetscReal,spbas_matrix*);
95b2f1ab58SBarry Smith PetscErrorCode spbas_matrix_to_crs(spbas_matrix, MatScalar **,PetscInt **,PetscInt**);
96b2f1ab58SBarry Smith PetscErrorCode spbas_dump(const char*,spbas_matrix);
97b2f1ab58SBarry Smith PetscErrorCode spbas_transpose(spbas_matrix,spbas_matrix*);
98b2f1ab58SBarry Smith PetscErrorCode spbas_apply_reordering(spbas_matrix*, const PetscInt *, const PetscInt*);
99b2f1ab58SBarry Smith PetscErrorCode spbas_pattern_only(PetscInt, PetscInt, PetscInt *, PetscInt*, spbas_matrix *);
100b2f1ab58SBarry Smith PetscErrorCode spbas_power (spbas_matrix, PetscInt, spbas_matrix *);
101b2f1ab58SBarry Smith PetscErrorCode spbas_keep_upper(spbas_matrix *);
102b2f1ab58SBarry Smith 
103