Actual source code: spbas.h

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: /*
  2:    Define type spbas_matrix: sparse matrices using pointers

  4:    Global matrix information
  5:       nrows, ncols: dimensions
  6:       nnz         : number of nonzeros (in whole matrix)
  7:       col_idx_type: storage scheme for column numbers
  8:                     SPBAS_COLUMN_NUMBERS:
  9:                         array icol contains column indices:
 10:                            A(i,icol[i][j]) = values[i][j]
 11:                     SPBAS_DIAGONAL_OFFSETS:
 12:                         array icol contains diagonal offsets:
 13:                            A(i,i+icol[i][j]) = values[i][j]
 14:                     SPBAS_OFFSET_ARRAY:
 15:                         array icol contains offsets wrt array
 16:                         icol0:
 17:                            A(i,icol0[i]+icol[i][j]) = values[i][j]

 19:    Information about each row
 20:       row_nnz     : number of nonzeros for each row
 21:       icol0       : column index offset (when needed, otherwise NULL)
 22:       icols       : array of diagonal offsets for each row, as descibed
 23:                     for col_idx_type, above
 24:       values      : array of matrix entries for each row
 25:                     when values == NULL, this matrix is really
 26:                     a sparseness pattern, not a matrix

 28:    The other fields describe the way in which the data are stored
 29:    in memory.

 31:       block_data  : The pointers icols[i] all point to places in a
 32:                     single allocated array. Only for icols[0] was
 33:                     malloc called. Freeing icols[0] will free
 34:                     all other icols=arrays as well.
 35:                     Same for arrays values[i]
 36: */

 38: #define SPBAS_COLUMN_NUMBERS   (0)
 39: #define SPBAS_DIAGONAL_OFFSETS (1)
 40: #define SPBAS_OFFSET_ARRAY     (2)

 42: #define NEGATIVE_DIAGONAL (-42)

 44: typedef struct {
 45:   PetscInt nrows;
 46:   PetscInt ncols;
 47:   PetscInt nnz;
 48:   PetscInt col_idx_type;

 50:   PetscInt    *row_nnz;
 51:   PetscInt    *icol0;
 52:   PetscInt    **icols;
 53:   PetscScalar **values;

 55:   PetscBool   block_data;
 56:   PetscInt    n_alloc_icol;
 57:   PetscInt    n_alloc_val;
 58:   PetscInt    *alloc_icol;
 59:   PetscScalar *alloc_val;
 60: } spbas_matrix;


 63: /*
 64:   spbas_compress_pattern:
 65:      calculate a compressed sparseness pattern for a sparseness pattern
 66:      given in compressed row storage. The compressed sparseness pattern may
 67:      require (much) less memory.

 69:   spbas_memory_requirement:
 70:      Calculate the number of bytes needed to store tha matrix

 72:   spbas_incomplete_cholesky:
 73:      Incomplete Cholesky decomposition

 75:   spbas_delete:
 76:      de-allocate the arrays owned by this matrix

 78:   spbas_matrix_to_crs:
 79:      Convert an spbas_matrix to compessed row storage

 81:   spbas_dump:
 82:      Print the matrix in i,j,val-format

 84:   spbas_transpose:
 85:      Return the transpose of a matrix

 87:   spbas_pattern_only:
 88:      Return the sparseness pattern (matrix without values) of a
 89:      compressed row storage
 90: */
 91: PetscErrorCode spbas_compress_pattern(PetscInt*,PetscInt*,PetscInt,PetscInt,PetscInt,spbas_matrix*,PetscReal*);
 92: long int       spbas_memory_requirement(spbas_matrix);
 93: PetscErrorCode spbas_delete(spbas_matrix);
 94: PetscErrorCode spbas_incomplete_cholesky(Mat,const PetscInt*,const PetscInt*,spbas_matrix,PetscReal,PetscReal,spbas_matrix*);
 95: PetscErrorCode spbas_matrix_to_crs(spbas_matrix, MatScalar **,PetscInt **,PetscInt**);
 96: PetscErrorCode spbas_dump(const char*,spbas_matrix);
 97: PetscErrorCode spbas_transpose(spbas_matrix,spbas_matrix*);
 98: PetscErrorCode spbas_apply_reordering(spbas_matrix*, const PetscInt*, const PetscInt*);
 99: PetscErrorCode spbas_pattern_only(PetscInt, PetscInt, PetscInt*, PetscInt*, spbas_matrix*);
100: PetscErrorCode spbas_power (spbas_matrix, PetscInt, spbas_matrix*);
101: PetscErrorCode spbas_keep_upper(spbas_matrix*);