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