Blame spmatrix/gsl_spmatrix.h

Packit 67cb25
/* gsl_spmatrix.h
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 2012-2014 Patrick Alken
Packit 67cb25
 * 
Packit 67cb25
 * This program is free software; you can redistribute it and/or modify
Packit 67cb25
 * it under the terms of the GNU General Public License as published by
Packit 67cb25
 * the Free Software Foundation; either version 3 of the License, or (at
Packit 67cb25
 * your option) any later version.
Packit 67cb25
 * 
Packit 67cb25
 * This program is distributed in the hope that it will be useful, but
Packit 67cb25
 * WITHOUT ANY WARRANTY; without even the implied warranty of
Packit 67cb25
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Packit 67cb25
 * General Public License for more details.
Packit 67cb25
 * 
Packit 67cb25
 * You should have received a copy of the GNU General Public License
Packit 67cb25
 * along with this program; if not, write to the Free Software
Packit 67cb25
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
#ifndef __GSL_SPMATRIX_H__
Packit 67cb25
#define __GSL_SPMATRIX_H__
Packit 67cb25
Packit 67cb25
#include <stdlib.h>
Packit 67cb25
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_vector.h>
Packit 67cb25
#include <gsl/gsl_matrix.h>
Packit 67cb25
Packit 67cb25
#undef __BEGIN_DECLS
Packit 67cb25
#undef __END_DECLS
Packit 67cb25
#ifdef __cplusplus
Packit 67cb25
# define __BEGIN_DECLS extern "C" {
Packit 67cb25
# define __END_DECLS }
Packit 67cb25
#else
Packit 67cb25
# define __BEGIN_DECLS /* empty */
Packit 67cb25
# define __END_DECLS /* empty */
Packit 67cb25
#endif
Packit 67cb25
Packit 67cb25
__BEGIN_DECLS
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
 * Binary tree data structure for storing sparse matrix elements
Packit 67cb25
 * in triplet format. This is used for efficiently detecting
Packit 67cb25
 * duplicates and element retrieval via gsl_spmatrix_get
Packit 67cb25
 */
Packit 67cb25
typedef struct
Packit 67cb25
{
Packit 67cb25
  void *tree;       /* tree structure */
Packit 67cb25
  void *node_array; /* preallocated array of tree nodes */
Packit 67cb25
  size_t n;         /* number of tree nodes in use (<= nzmax) */
Packit 67cb25
} gsl_spmatrix_tree;
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
 * Triplet format:
Packit 67cb25
 *
Packit 67cb25
 * If data[n] = A_{ij}, then:
Packit 67cb25
 *   i = A->i[n]
Packit 67cb25
 *   j = A->p[n]
Packit 67cb25
 *
Packit 67cb25
 * Compressed column format (CCS):
Packit 67cb25
 *
Packit 67cb25
 * If data[n] = A_{ij}, then:
Packit 67cb25
 *   i = A->i[n]
Packit 67cb25
 *   A->p[j] <= n < A->p[j+1]
Packit 67cb25
 * so that column j is stored in
Packit 67cb25
 * [ data[p[j]], data[p[j] + 1], ..., data[p[j+1] - 1] ]
Packit 67cb25
 *
Packit 67cb25
 * Compressed row format (CRS):
Packit 67cb25
 *
Packit 67cb25
 * If data[n] = A_{ij}, then:
Packit 67cb25
 *   j = A->i[n]
Packit 67cb25
 *   A->p[i] <= n < A->p[i+1]
Packit 67cb25
 * so that row i is stored in
Packit 67cb25
 * [ data[p[i]], data[p[i] + 1], ..., data[p[i+1] - 1] ]
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
typedef struct
Packit 67cb25
{
Packit 67cb25
  size_t size1;  /* number of rows */
Packit 67cb25
  size_t size2;  /* number of columns */
Packit 67cb25
Packit 67cb25
  /* i (size nzmax) contains:
Packit 67cb25
   *
Packit 67cb25
   * Triplet/CCS: row indices
Packit 67cb25
   * CRS: column indices
Packit 67cb25
   */
Packit 67cb25
  size_t *i;
Packit 67cb25
Packit 67cb25
  double *data;  /* matrix elements of size nzmax */
Packit 67cb25
Packit 67cb25
  /*
Packit 67cb25
   * p contains the column indices (triplet) or column pointers (compcol)
Packit 67cb25
   *
Packit 67cb25
   * triplet: p[n] = column number of element data[n]
Packit 67cb25
   * CCS:     p[j] = index in data of first non-zero element in column j
Packit 67cb25
   * CRS:     p[i] = index in data of first non-zero element in row i
Packit 67cb25
   */
Packit 67cb25
  size_t *p;
Packit 67cb25
Packit 67cb25
  size_t nzmax;  /* maximum number of matrix elements */
Packit 67cb25
  size_t nz;     /* number of non-zero values in matrix */
Packit 67cb25
Packit 67cb25
  gsl_spmatrix_tree *tree_data; /* binary tree for sorting triplet data */
Packit 67cb25
Packit 67cb25
  /*
Packit 67cb25
   * workspace of size MAX(size1,size2)*MAX(sizeof(double),sizeof(size_t))
Packit 67cb25
   * used in various routines
Packit 67cb25
   */
Packit 67cb25
  union
Packit 67cb25
    {
Packit 67cb25
      void *work;
Packit 67cb25
      size_t *work_sze;
Packit 67cb25
      double *work_dbl;
Packit 67cb25
    };
Packit 67cb25
Packit 67cb25
  size_t sptype; /* sparse storage type */
Packit 67cb25
} gsl_spmatrix;
Packit 67cb25
Packit 67cb25
#define GSL_SPMATRIX_TRIPLET      (0)
Packit 67cb25
#define GSL_SPMATRIX_CCS          (1)
Packit 67cb25
#define GSL_SPMATRIX_CRS          (2)
Packit 67cb25
Packit 67cb25
#define GSL_SPMATRIX_ISTRIPLET(m) ((m)->sptype == GSL_SPMATRIX_TRIPLET)
Packit 67cb25
#define GSL_SPMATRIX_ISCCS(m)     ((m)->sptype == GSL_SPMATRIX_CCS)
Packit 67cb25
#define GSL_SPMATRIX_ISCRS(m)     ((m)->sptype == GSL_SPMATRIX_CRS)
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
 * Prototypes
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
gsl_spmatrix *gsl_spmatrix_alloc(const size_t n1, const size_t n2);
Packit 67cb25
gsl_spmatrix *gsl_spmatrix_alloc_nzmax(const size_t n1, const size_t n2,
Packit 67cb25
                                       const size_t nzmax, const size_t flags);
Packit 67cb25
void gsl_spmatrix_free(gsl_spmatrix *m);
Packit 67cb25
int gsl_spmatrix_realloc(const size_t nzmax, gsl_spmatrix *m);
Packit 67cb25
int gsl_spmatrix_set_zero(gsl_spmatrix *m);
Packit 67cb25
size_t gsl_spmatrix_nnz(const gsl_spmatrix *m);
Packit 67cb25
int gsl_spmatrix_compare_idx(const size_t ia, const size_t ja,
Packit 67cb25
                             const size_t ib, const size_t jb);
Packit 67cb25
int gsl_spmatrix_tree_rebuild(gsl_spmatrix * m);
Packit 67cb25
Packit 67cb25
/* spcopy.c */
Packit 67cb25
int gsl_spmatrix_memcpy(gsl_spmatrix *dest, const gsl_spmatrix *src);
Packit 67cb25
Packit 67cb25
/* spgetset.c */
Packit 67cb25
double gsl_spmatrix_get(const gsl_spmatrix *m, const size_t i,
Packit 67cb25
                        const size_t j);
Packit 67cb25
int gsl_spmatrix_set(gsl_spmatrix *m, const size_t i, const size_t j,
Packit 67cb25
                     const double x);
Packit 67cb25
double *gsl_spmatrix_ptr(gsl_spmatrix *m, const size_t i, const size_t j);
Packit 67cb25
Packit 67cb25
/* spcompress.c */
Packit 67cb25
gsl_spmatrix *gsl_spmatrix_compcol(const gsl_spmatrix *T);
Packit 67cb25
gsl_spmatrix *gsl_spmatrix_ccs(const gsl_spmatrix *T);
Packit 67cb25
gsl_spmatrix *gsl_spmatrix_crs(const gsl_spmatrix *T);
Packit 67cb25
void gsl_spmatrix_cumsum(const size_t n, size_t *c);
Packit 67cb25
Packit 67cb25
/* spio.c */
Packit 67cb25
int gsl_spmatrix_fprintf(FILE *stream, const gsl_spmatrix *m,
Packit 67cb25
                         const char *format);
Packit 67cb25
gsl_spmatrix * gsl_spmatrix_fscanf(FILE *stream);
Packit 67cb25
int gsl_spmatrix_fwrite(FILE *stream, const gsl_spmatrix *m);
Packit 67cb25
int gsl_spmatrix_fread(FILE *stream, gsl_spmatrix *m);
Packit 67cb25
Packit 67cb25
/* spoper.c */
Packit 67cb25
int gsl_spmatrix_scale(gsl_spmatrix *m, const double x);
Packit 67cb25
int gsl_spmatrix_minmax(const gsl_spmatrix *m, double *min_out,
Packit 67cb25
                        double *max_out);
Packit 67cb25
int gsl_spmatrix_add(gsl_spmatrix *c, const gsl_spmatrix *a,
Packit 67cb25
                     const gsl_spmatrix *b);
Packit 67cb25
int gsl_spmatrix_d2sp(gsl_spmatrix *S, const gsl_matrix *A);
Packit 67cb25
int gsl_spmatrix_sp2d(gsl_matrix *A, const gsl_spmatrix *S);
Packit 67cb25
Packit 67cb25
/* spprop.c */
Packit 67cb25
int gsl_spmatrix_equal(const gsl_spmatrix *a, const gsl_spmatrix *b);
Packit 67cb25
Packit 67cb25
/* spswap.c */
Packit 67cb25
int gsl_spmatrix_transpose(gsl_spmatrix * m);
Packit 67cb25
int gsl_spmatrix_transpose2(gsl_spmatrix * m);
Packit 67cb25
int gsl_spmatrix_transpose_memcpy(gsl_spmatrix *dest, const gsl_spmatrix *src);
Packit 67cb25
Packit 67cb25
__END_DECLS
Packit 67cb25
Packit 67cb25
#endif /* __GSL_SPMATRIX_H__ */