cs.h
Go to the documentation of this file.
1 #ifndef _CS_H
2 #define _CS_H
3 #include <stdlib.h>
4 #include <math.h>
5 #include <stdio.h>
6 
7 #define CS_VER 2 /* CSparse Version 2.2.3 */
8 #define CS_SUBVER 2
9 #define CS_SUBSUB 3
10 #define CS_DATE "Jan 20, 2009" /* CSparse release date */
11 #define CS_COPYRIGHT "Copyright (c) Timothy A. Davis, 2006-2009"
12 
13 
14 
15 /* --- primary CSparse routines and data structures ------------------------- */
16 typedef struct cs_sparse /* matrix in compressed-column or triplet form */
17 {
18  int nzmax ; /* maximum number of entries */
19  int m ; /* number of rows */
20  int n ; /* number of columns */
21  int *p ; /* column pointers (size n+1) or col indices (size nzmax) */
22  int *i ; /* row indices, size nzmax */
23  double *x ; /* numerical values, size nzmax */
24  int nz ; /* # of entries in triplet matrix, -1 for compressed-col */
25 } cs ;
26 
27 cs *cs_add (const cs *A, const cs *B, double alpha, double beta) ;
28 int cs_cholsol (int order, const cs *A, double *b) ;
29 cs *cs_compress (const cs *T) ;
30 int cs_dupl (cs *A) ;
31 int cs_entry (cs *T, int i, int j, double x) ;
32 int cs_gaxpy (const cs *A, const double *x, double *y) ;
33 cs *cs_load (FILE *f) ;
34 int cs_lusol (int order, const cs *A, double *b, double tol) ;
35 cs *cs_multiply (const cs *A, const cs *B) ;
36 double cs_norm (const cs *A) ;
37 int cs_print (const cs *A, int brief) ;
38 int cs_qrsol (int order, const cs *A, double *b) ;
39 cs *cs_transpose (const cs *A, int values) ;
40 /* utilities */
41 void *cs_calloc (int n, size_t size) ;
42 void *cs_free (void *p) ;
43 void *cs_realloc (void *p, int n, size_t size, int *ok) ;
44 cs *cs_spalloc (int m, int n, int nzmax, int values, int triplet) ;
45 cs *cs_spfree (cs *A) ;
46 int cs_sprealloc (cs *A, int nzmax) ;
47 void *cs_malloc (int n, size_t size) ;
48 
49 /* --- secondary CSparse routines and data structures ----------------------- */
50 typedef struct cs_symbolic /* symbolic Cholesky, LU, or QR analysis */
51 {
52  int *pinv ; /* inverse row perm. for QR, fill red. perm for Chol */
53  int *q ; /* fill-reducing column permutation for LU and QR */
54  int *parent ; /* elimination tree for Cholesky and QR */
55  int *cp ; /* column pointers for Cholesky, row counts for QR */
56  int *leftmost ; /* leftmost[i] = min(find(A(i,:))), for QR */
57  int m2 ; /* # of rows for QR, after adding fictitious rows */
58  double lnz ; /* # entries in L for LU or Cholesky; in V for QR */
59  double unz ; /* # entries in U for LU; in R for QR */
60 } css ;
61 
62 typedef struct cs_numeric /* numeric Cholesky, LU, or QR factorization */
63 {
64  cs *L ; /* L for LU and Cholesky, V for QR */
65  cs *U ; /* U for LU, R for QR, not used for Cholesky */
66  int *pinv ; /* partial pivoting for LU */
67  double *B ; /* beta [0..n-1] for QR */
68 } csn ;
69 
70 typedef struct cs_dmperm_results /* cs_dmperm or cs_scc output */
71 {
72  int *p ; /* size m, row permutation */
73  int *q ; /* size n, column permutation */
74  int *r ; /* size nb+1, block k is rows r[k] to r[k+1]-1 in A(p,q) */
75  int *s ; /* size nb+1, block k is cols s[k] to s[k+1]-1 in A(p,q) */
76  int nb ; /* # of blocks in fine dmperm decomposition */
77  int rr [5] ; /* coarse row decomposition */
78  int cc [5] ; /* coarse column decomposition */
79 } csd ;
80 
81 int *cs_amd (int order, const cs *A) ;
82 csn *cs_chol (const cs *A, const css *S) ;
83 csd *cs_dmperm (const cs *A, int seed) ;
84 int cs_droptol (cs *A, double tol) ;
85 int cs_dropzeros (cs *A) ;
86 int cs_happly (const cs *V, int i, double beta, double *x) ;
87 int cs_ipvec (const int *p, const double *b, double *x, int n) ;
88 int cs_lsolve (const cs *L, double *x) ;
89 int cs_ltsolve (const cs *L, double *x) ;
90 csn *cs_lu (const cs *A, const css *S, double tol) ;
91 cs *cs_permute (const cs *A, const int *pinv, const int *q, int values) ;
92 int *cs_pinv (const int *p, int n) ;
93 int cs_pvec (const int *p, const double *b, double *x, int n) ;
94 csn *cs_qr (const cs *A, const css *S) ;
95 css *cs_schol (int order, const cs *A) ;
96 css *cs_sqr (int order, const cs *A, int qr) ;
97 cs *cs_symperm (const cs *A, const int *pinv, int values) ;
98 int cs_updown (cs *L, int sigma, const cs *C, const int *parent) ;
99 int cs_usolve (const cs *U, double *x) ;
100 int cs_utsolve (const cs *U, double *x) ;
101 /* utilities */
102 css *cs_sfree (css *S) ;
103 csn *cs_nfree (csn *N) ;
104 csd *cs_dfree (csd *D) ;
105 
106 /* --- tertiary CSparse routines -------------------------------------------- */
107 int *cs_counts (const cs *A, const int *parent, const int *post, int ata) ;
108 double cs_cumsum (int *p, int *c, int n) ;
109 int cs_dfs (int j, cs *G, int top, int *xi, int *pstack, const int *pinv) ;
110 int cs_ereach (const cs *A, int k, const int *parent, int *s, int *w) ;
111 int *cs_etree (const cs *A, int ata) ;
112 int cs_fkeep (cs *A, int (*fkeep) (int, int, double, void *), void *other) ;
113 double cs_house (double *x, double *beta, int n) ;
114 int cs_leaf (int i, int j, const int *first, int *maxfirst, int *prevleaf,
115  int *ancestor, int *jleaf) ;
116 int *cs_maxtrans (const cs *A, int seed) ;
117 int *cs_post (const int *parent, int n) ;
118 int *cs_randperm (int n, int seed) ;
119 int cs_reach (cs *G, const cs *B, int k, int *xi, const int *pinv) ;
120 int cs_scatter (const cs *A, int j, double beta, int *w, double *x, int mark,
121  cs *C, int nz) ;
122 csd *cs_scc (cs *A) ;
123 int cs_spsolve (cs *G, const cs *B, int k, int *xi, double *x,
124  const int *pinv, int lo) ;
125 int cs_tdfs (int j, int k, int *head, const int *next, int *post,
126  int *stack) ;
127 /* utilities */
128 csd *cs_dalloc (int m, int n) ;
129 csd *cs_ddone (csd *D, cs *C, void *w, int ok) ;
130 cs *cs_done (cs *C, void *w, void *x, int ok) ;
131 int *cs_idone (int *p, cs *C, void *w, int ok) ;
132 csn *cs_ndone (csn *N, cs *C, void *w, void *x, int ok) ;
133 
134 #define CS_MAX(a,b) (((a) > (b)) ? (a) : (b))
135 #define CS_MIN(a,b) (((a) < (b)) ? (a) : (b))
136 #define CS_FLIP(i) (-(i)-2)
137 #define CS_UNFLIP(i) (((i) < 0) ? CS_FLIP(i) : (i))
138 #define CS_MARKED(w,j) (w [j] < 0)
139 #define CS_MARK(w,j) { w [j] = CS_FLIP (w [j]) ; }
140 #define CS_CSC(A) (A && (A->nz == -1))
141 #define CS_TRIPLET(A) (A && (A->nz >= 0))
142 #endif
#define N
css * cs_sqr(int order, const cs *A, int qr)
Definition: cs_sqr.c:60
int cs_lusol(int order, const cs *A, double *b, double tol)
Definition: cs_lusol.c:3
int cs_reach(cs *G, const cs *B, int k, int *xi, const int *pinv)
Definition: cs_reach.c:4
int cs_dupl(cs *A)
Definition: cs_dupl.c:3
cs * cs_load(FILE *f)
Definition: cs_load.c:3
USING_NAMESPACE_ACADO typedef TaylorVariable< Interval > T
csn * cs_chol(const cs *A, const css *S)
Definition: cs_chol.c:3
int * pinv
Definition: cs.h:66
int cs_lsolve(const cs *L, double *x)
Definition: cs_lsolve.c:3
struct cs_dmperm_results csd
double * B
Definition: cs.h:67
csn * cs_nfree(csn *N)
Definition: cs_util.c:42
cs * cs_spfree(cs *A)
Definition: cs_util.c:32
cs * cs_spalloc(int m, int n, int nzmax, int values, int triplet)
Definition: cs_util.c:3
int cs_ltsolve(const cs *L, double *x)
Definition: cs_ltsolve.c:3
csn * cs_qr(const cs *A, const css *S)
Definition: cs_qr.c:3
int * s
Definition: cs.h:75
cs * cs_done(cs *C, void *w, void *x, int ok)
Definition: cs_util.c:89
cs * cs_symperm(const cs *A, const int *pinv, int values)
Definition: cs_symperm.c:3
cs * cs_permute(const cs *A, const int *pinv, const int *q, int values)
Definition: cs_permute.c:3
int cs_droptol(cs *A, double tol)
Definition: cs_droptol.c:6
int n
Definition: cs.h:20
cs * U
Definition: cs.h:65
csd * cs_dalloc(int m, int n)
Definition: cs_util.c:65
csd * cs_dmperm(const cs *A, int seed)
Definition: cs_dmperm.c:68
void * cs_free(void *p)
Definition: cs_malloc.c:22
int * cs_post(const int *parent, int n)
Definition: cs_post.c:3
int cs_fkeep(cs *A, int(*fkeep)(int, int, double, void *), void *other)
Definition: cs_fkeep.c:3
struct cs_numeric csn
double cs_house(double *x, double *beta, int n)
Definition: cs_house.c:4
double lnz
Definition: cs.h:58
int cs_qrsol(int order, const cs *A, double *b)
Definition: cs_qrsol.c:3
Definition: cs.h:62
int cs_ipvec(const int *p, const double *b, double *x, int n)
Definition: cs_ipvec.c:3
cs * cs_transpose(const cs *A, int values)
Definition: cs_transpose.c:3
int nzmax
Definition: cs.h:18
int cs_usolve(const cs *U, double *x)
Definition: cs_usolve.c:3
int cs_dfs(int j, cs *G, int top, int *xi, int *pstack, const int *pinv)
Definition: cs_dfs.c:3
int cs_print(const cs *A, int brief)
Definition: cs_print.c:3
int cs_happly(const cs *V, int i, double beta, double *x)
Definition: cs_happly.c:3
int cs_sprealloc(cs *A, int nzmax)
Definition: cs_util.c:18
int cs_ereach(const cs *A, int k, const int *parent, int *s, int *w)
Definition: cs_ereach.c:3
double unz
Definition: cs.h:59
csd * cs_dfree(csd *D)
Definition: cs_util.c:78
int * cs_etree(const cs *A, int ata)
Definition: cs_etree.c:3
struct cs_sparse cs
int nz
Definition: cs.h:24
int cs_utsolve(const cs *U, double *x)
Definition: cs_utsolve.c:3
int cs_leaf(int i, int j, const int *first, int *maxfirst, int *prevleaf, int *ancestor, int *jleaf)
Definition: cs_leaf.c:3
struct cs_symbolic css
int * cs_maxtrans(const cs *A, int seed)
Definition: cs_maxtrans.c:44
int * leftmost
Definition: cs.h:56
int * cs_idone(int *p, cs *C, void *w, int ok)
Definition: cs_util.c:97
int cs_entry(cs *T, int i, int j, double x)
Definition: cs_entry.c:3
cs * cs_add(const cs *A, const cs *B, double alpha, double beta)
Definition: cs_add.c:3
int * p
Definition: cs.h:21
double cs_norm(const cs *A)
Definition: cs_norm.c:3
int cs_scatter(const cs *A, int j, double beta, int *w, double *x, int mark, cs *C, int nz)
Definition: cs_scatter.c:3
Definition: cs.h:16
int * i
Definition: cs.h:22
int * cp
Definition: cs.h:55
void * cs_calloc(int n, size_t size)
Definition: cs_malloc.c:16
int cs_pvec(const int *p, const double *b, double *x, int n)
Definition: cs_pvec.c:3
int * pinv
Definition: cs.h:52
double cs_cumsum(int *p, int *c, int n)
Definition: cs_cumsum.c:3
cs * cs_multiply(const cs *A, const cs *B)
Definition: cs_multiply.c:3
int * q
Definition: cs.h:53
int cs_tdfs(int j, int k, int *head, const int *next, int *post, int *stack)
Definition: cs_tdfs.c:3
cs * cs_compress(const cs *T)
Definition: cs_compress.c:3
csd * cs_ddone(csd *D, cs *C, void *w, int ok)
Definition: cs_util.c:114
void * cs_malloc(int n, size_t size)
Definition: cs_malloc.c:10
csd * cs_scc(cs *A)
Definition: cs_scc.c:3
#define L
SegmentReturnType head(Index vecSize)
Definition: BlockMethods.h:781
int cs_cholsol(int order, const cs *A, double *b)
Definition: cs_cholsol.c:3
void * cs_realloc(void *p, int n, size_t size, int *ok)
Definition: cs_malloc.c:29
int * cs_amd(int order, const cs *A)
Definition: cs_amd.c:18
int * p
Definition: cs.h:72
int * r
Definition: cs.h:74
int cs_dropzeros(cs *A)
Definition: cs_dropzeros.c:6
int m2
Definition: cs.h:57
Expression next(const Expression &arg)
int * cs_pinv(const int *p, int n)
Definition: cs_pinv.c:3
double * x
Definition: cs.h:23
css * cs_sfree(css *S)
Definition: cs_util.c:53
int * parent
Definition: cs.h:54
csn * cs_ndone(csn *N, cs *C, void *w, void *x, int ok)
Definition: cs_util.c:105
int cs_gaxpy(const cs *A, const double *x, double *y)
Definition: cs_gaxpy.c:3
csn * cs_lu(const cs *A, const css *S, double tol)
Definition: cs_lu.c:3
css * cs_schol(int order, const cs *A)
Definition: cs_schol.c:3
int cs_spsolve(cs *G, const cs *B, int k, int *xi, double *x, const int *pinv, int lo)
Definition: cs_spsolve.c:3
int * cs_randperm(int n, int seed)
Definition: cs_randperm.c:5
int cs_updown(cs *L, int sigma, const cs *C, const int *parent)
Definition: cs_updown.c:3
int * q
Definition: cs.h:73
int m
Definition: cs.h:19
int * cs_counts(const cs *A, const int *parent, const int *post, int ata)
Definition: cs_counts.c:17
cs * L
Definition: cs.h:64
Definition: cs.h:50


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:34:31