$search
00001 #ifndef _CS_H 00002 #define _CS_H 00003 #include <stdlib.h> 00004 #include <limits.h> 00005 #include <math.h> 00006 #include <stdio.h> 00007 #ifdef MATLAB_MEX_FILE 00008 #include "mex.h" 00009 #endif 00010 #define CS_VER 2 /* CSparse Version 2.2.3 */ 00011 #define CS_SUBVER 2 00012 #define CS_SUBSUB 3 00013 #define CS_DATE "Jan 20, 2009" /* CSparse release date */ 00014 #define CS_COPYRIGHT "Copyright (c) Timothy A. Davis, 2006-2009" 00015 00016 /* --- primary CSparse routines and data structures ------------------------- */ 00017 typedef struct cs_sparse /* matrix in compressed-column or triplet form */ 00018 { 00019 int nzmax ; /* maximum number of entries */ 00020 int m ; /* number of rows */ 00021 int n ; /* number of columns */ 00022 int *p ; /* column pointers (size n+1) or col indices (size nzmax) */ 00023 int *i ; /* row indices, size nzmax */ 00024 double *x ; /* numerical values, size nzmax */ 00025 int nz ; /* # of entries in triplet matrix, -1 for compressed-col */ 00026 } cs ; 00027 00028 cs *cs_add (const cs *A, const cs *B, double alpha, double beta) ; 00029 int cs_cholsol (int order, const cs *A, double *b) ; 00030 cs *cs_compress (const cs *T) ; 00031 int cs_dupl (cs *A) ; 00032 int cs_entry (cs *T, int i, int j, double x) ; 00033 int cs_gaxpy (const cs *A, const double *x, double *y) ; 00034 cs *cs_load (FILE *f) ; 00035 int cs_lusol (int order, const cs *A, double *b, double tol) ; 00036 cs *cs_multiply (const cs *A, const cs *B) ; 00037 double cs_norm (const cs *A) ; 00038 int cs_print (const cs *A, int brief) ; 00039 int cs_qrsol (int order, const cs *A, double *b) ; 00040 cs *cs_transpose (const cs *A, int values) ; 00041 /* utilities */ 00042 void *cs_calloc (int n, size_t size) ; 00043 void *cs_free (void *p) ; 00044 void *cs_realloc (void *p, int n, size_t size, int *ok) ; 00045 cs *cs_spalloc (int m, int n, int nzmax, int values, int triplet) ; 00046 cs *cs_spfree (cs *A) ; 00047 int cs_sprealloc (cs *A, int nzmax) ; 00048 void *cs_malloc (int n, size_t size) ; 00049 00050 /* --- secondary CSparse routines and data structures ----------------------- */ 00051 typedef struct cs_symbolic /* symbolic Cholesky, LU, or QR analysis */ 00052 { 00053 int *pinv ; /* inverse row perm. for QR, fill red. perm for Chol */ 00054 int *q ; /* fill-reducing column permutation for LU and QR */ 00055 int *parent ; /* elimination tree for Cholesky and QR */ 00056 int *cp ; /* column pointers for Cholesky, row counts for QR */ 00057 int *leftmost ; /* leftmost[i] = min(find(A(i,:))), for QR */ 00058 int m2 ; /* # of rows for QR, after adding fictitious rows */ 00059 double lnz ; /* # entries in L for LU or Cholesky; in V for QR */ 00060 double unz ; /* # entries in U for LU; in R for QR */ 00061 } css ; 00062 00063 typedef struct cs_numeric /* numeric Cholesky, LU, or QR factorization */ 00064 { 00065 cs *L ; /* L for LU and Cholesky, V for QR */ 00066 cs *U ; /* U for LU, R for QR, not used for Cholesky */ 00067 int *pinv ; /* partial pivoting for LU */ 00068 double *B ; /* beta [0..n-1] for QR */ 00069 } csn ; 00070 00071 typedef struct cs_dmperm_results /* cs_dmperm or cs_scc output */ 00072 { 00073 int *p ; /* size m, row permutation */ 00074 int *q ; /* size n, column permutation */ 00075 int *r ; /* size nb+1, block k is rows r[k] to r[k+1]-1 in A(p,q) */ 00076 int *s ; /* size nb+1, block k is cols s[k] to s[k+1]-1 in A(p,q) */ 00077 int nb ; /* # of blocks in fine dmperm decomposition */ 00078 int rr [5] ; /* coarse row decomposition */ 00079 int cc [5] ; /* coarse column decomposition */ 00080 } csd ; 00081 00082 int *cs_amd (int order, const cs *A) ; 00083 csn *cs_chol (const cs *A, const css *S) ; 00084 csd *cs_dmperm (const cs *A, int seed) ; 00085 int cs_droptol (cs *A, double tol) ; 00086 int cs_dropzeros (cs *A) ; 00087 int cs_happly (const cs *V, int i, double beta, double *x) ; 00088 int cs_ipvec (const int *p, const double *b, double *x, int n) ; 00089 int cs_lsolve (const cs *L, double *x) ; 00090 int cs_ltsolve (const cs *L, double *x) ; 00091 csn *cs_lu (const cs *A, const css *S, double tol) ; 00092 cs *cs_permute (const cs *A, const int *pinv, const int *q, int values) ; 00093 int *cs_pinv (const int *p, int n) ; 00094 int cs_pvec (const int *p, const double *b, double *x, int n) ; 00095 csn *cs_qr (const cs *A, const css *S) ; 00096 css *cs_schol (int order, const cs *A) ; 00097 css *cs_sqr (int order, const cs *A, int qr) ; 00098 cs *cs_symperm (const cs *A, const int *pinv, int values) ; 00099 int cs_updown (cs *L, int sigma, const cs *C, const int *parent) ; 00100 int cs_usolve (const cs *U, double *x) ; 00101 int cs_utsolve (const cs *U, double *x) ; 00102 /* utilities */ 00103 css *cs_sfree (css *S) ; 00104 csn *cs_nfree (csn *N) ; 00105 csd *cs_dfree (csd *D) ; 00106 00107 /* --- tertiary CSparse routines -------------------------------------------- */ 00108 int *cs_counts (const cs *A, const int *parent, const int *post, int ata) ; 00109 double cs_cumsum (int *p, int *c, int n) ; 00110 int cs_dfs (int j, cs *G, int top, int *xi, int *pstack, const int *pinv) ; 00111 int cs_ereach (const cs *A, int k, const int *parent, int *s, int *w) ; 00112 int *cs_etree (const cs *A, int ata) ; 00113 int cs_fkeep (cs *A, int (*fkeep) (int, int, double, void *), void *other) ; 00114 double cs_house (double *x, double *beta, int n) ; 00115 int cs_leaf (int i, int j, const int *first, int *maxfirst, int *prevleaf, 00116 int *ancestor, int *jleaf) ; 00117 int *cs_maxtrans (const cs *A, int seed) ; 00118 int *cs_post (const int *parent, int n) ; 00119 int *cs_randperm (int n, int seed) ; 00120 int cs_reach (cs *G, const cs *B, int k, int *xi, const int *pinv) ; 00121 int cs_scatter (const cs *A, int j, double beta, int *w, double *x, int mark, 00122 cs *C, int nz) ; 00123 csd *cs_scc (cs *A) ; 00124 int cs_spsolve (cs *G, const cs *B, int k, int *xi, double *x, 00125 const int *pinv, int lo) ; 00126 int cs_tdfs (int j, int k, int *head, const int *next, int *post, 00127 int *stack) ; 00128 /* utilities */ 00129 csd *cs_dalloc (int m, int n) ; 00130 csd *cs_ddone (csd *D, cs *C, void *w, int ok) ; 00131 cs *cs_done (cs *C, void *w, void *x, int ok) ; 00132 int *cs_idone (int *p, cs *C, void *w, int ok) ; 00133 csn *cs_ndone (csn *N, cs *C, void *w, void *x, int ok) ; 00134 00135 #define CS_MAX(a,b) (((a) > (b)) ? (a) : (b)) 00136 #define CS_MIN(a,b) (((a) < (b)) ? (a) : (b)) 00137 #define CS_FLIP(i) (-(i)-2) 00138 #define CS_UNFLIP(i) (((i) < 0) ? CS_FLIP(i) : (i)) 00139 #define CS_MARKED(w,j) (w [j] < 0) 00140 #define CS_MARK(w,j) { w [j] = CS_FLIP (w [j]) ; } 00141 #define CS_CSC(A) (A && (A->nz == -1)) 00142 #define CS_TRIPLET(A) (A && (A->nz >= 0)) 00143 00144 #endif