egsl_ops.c
Go to the documentation of this file.
1 #include <gsl/gsl_matrix.h>
2 #include <gsl/gsl_blas.h>
3 #include <gsl/gsl_linalg.h>
4 #include <gsl/gsl_eigen.h>
5 
6 #include "egsl.h"
7 
8 
9 
11  return egsl_sum(v1, egsl_scale(-1.0,v2));
12 }
13 
15  gsl_matrix *m1 = egsl_gslm(v1);
16  gsl_matrix *m2 = egsl_gslm(v2);
17  egsl_expect_size(v2, 0, m1->size2);
18  val v3 = egsl_alloc(m1->size1+m2->size1,m1->size2);
19  gsl_matrix *m3 = egsl_gslm(v3);
20  size_t i,j;
21  for(j=0;j<m1->size2;j++) {
22  for(i=0;i<m1->size1;i++)
23  gsl_matrix_set(m3, i, j, gsl_matrix_get(m1,i,j));
24 
25  for(i=0;i<m2->size1;i++)
26  gsl_matrix_set(m3, m1->size1+i, j, gsl_matrix_get(m2,i,j));
27  }
28  return v3;
29 }
30 
32  gsl_matrix *m1 = egsl_gslm(v1);
33  gsl_matrix *m2 = egsl_gslm(v2);
34  egsl_expect_size(v2, m1->size1, 0);
35  val v3 = egsl_alloc(m1->size1, m1->size2 + m2->size2);
36  gsl_matrix *m3 = egsl_gslm(v3);
37  size_t i,j;
38  for(i=0;i<m1->size1;i++) {
39  for(j=0;j<m1->size2;j++)
40  gsl_matrix_set(m3, i, j, gsl_matrix_get(m1,i,j));
41 
42  for(j=0;j<m2->size2;j++)
43  gsl_matrix_set(m3, i, m1->size2+j, gsl_matrix_get(m2,i,j));
44  }
45  return v3;
46 }
47 
48 void egsl_add_to(val v1, val v2) {
49  gsl_matrix * m1 = egsl_gslm(v1);
50  gsl_matrix * m2 = egsl_gslm(v2);
51  gsl_matrix_add(m1,m2);
52 }
53 
54 void egsl_add_to_col(val v1, size_t j, val v2) {
55 /* egsl_print("m1",v1);
56  egsl_print("m2",v2); */
57  gsl_matrix * m1 = egsl_gslm(v1);
58  gsl_matrix * m2 = egsl_gslm(v2);
59 
60 /* printf("m1 size = %d,%d j = %d\n",m1->size1,m1->size2,j); */
61  egsl_expect_size(v2, m1->size1, 1);
62  size_t i;
63  for(i=0;i<m1->size1;i++) {
64  *gsl_matrix_ptr(m1, i, j) += gsl_matrix_get(m2,i,0);
65  }
66 }
67 
68 
70  gsl_matrix * m1 = egsl_gslm(v1);
71  val v2 = egsl_alloc(m1->size1,m1->size2);
72  gsl_matrix * m2 = egsl_gslm(v2);
73  gsl_matrix_memcpy(m2,m1);
74  return v2;
75 }
76 
77 val egsl_scale(double s, val v1){
78  val v2 = egsl_copy_val(v1);
79  gsl_matrix * m2 = egsl_gslm(v2);
80  gsl_matrix_scale(m2, s);
81  return v2;
82 }
83 
84 val egsl_sum(val v1, val v2){
85  gsl_matrix * m1 = egsl_gslm(v1);
86  gsl_matrix * m2 = egsl_gslm(v2);
87  val v3 = egsl_alloc(m1->size1,m1->size2);
88  gsl_matrix * m3 = egsl_gslm(v3);
89  gsl_matrix_memcpy(m3,m1);
90  gsl_matrix_add(m3,m2);
91  return v3;
92 }
93 
94 val egsl_sum3(val v1, val v2, val v3){
95  return egsl_sum(v1, egsl_sum(v2,v3));
96 }
97 
98 val egsl_mult(val v1, val v2){
99  gsl_matrix * a = egsl_gslm(v1);
100  gsl_matrix * b = egsl_gslm(v2);
101  val v = egsl_alloc(a->size1,b->size2);
102  gsl_matrix * ab = egsl_gslm(v);
103  gsl_blas_dgemm(CblasNoTrans,CblasNoTrans,1.0,a,b,0.0,ab);
104  return v;
105 }
106 
108  gsl_matrix * m1 = egsl_gslm(v1);
109  val v2 = egsl_alloc(m1->size2,m1->size1);
110  gsl_matrix * m2 = egsl_gslm(v2);
111  gsl_matrix_transpose_memcpy(m2,m1);
112  return v2;
113 }
114 
116  gsl_matrix*A = egsl_gslm(v1);
117  val v2 = egsl_alloc(A->size1,A->size1);
118  gsl_matrix*invA = egsl_gslm(v2);
119  size_t n = A->size1;
120  gsl_matrix * m = gsl_matrix_alloc(n,n);
121  gsl_matrix_memcpy(m,A);
122  gsl_permutation * perm = gsl_permutation_alloc (n);
123  /* Make LU decomposition of matrix m */
124  int s;
125  gsl_linalg_LU_decomp (m, perm, &s);
126  /* Invert the matrix m */
127  gsl_linalg_LU_invert (m, perm, invA);
128  gsl_permutation_free(perm);
129  gsl_matrix_free(m);
130  return v2;
131 }
132 
133 void egsl_symm_eig(val v, double* eigenvalues, val* eigenvectors) {
134  gsl_matrix *m = egsl_gslm(v);
135  size_t N = m->size1;
136  /* Check for v to be square */
137 
138  gsl_matrix *A = gsl_matrix_alloc(N,N);
139  gsl_matrix_memcpy(A, m);
140 
141  gsl_vector *eval = gsl_vector_alloc(N);
142  gsl_matrix *evec = gsl_matrix_alloc(N,N);
143 
144  gsl_eigen_symmv_workspace * ws = gsl_eigen_symmv_alloc(N);
145  gsl_eigen_symmv(A, eval, evec, ws);
146  gsl_eigen_symmv_free(ws);
147 
148 
149  gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_VAL_DESC);
150 
151  size_t j;
152  for(j=0;j<N;j++) {
153  eigenvalues[j] = gsl_vector_get(eval, j);
154  eigenvectors[j] = egsl_alloc(N,1);
155  size_t i;
156  for(i=0;i<N;i++)
157  *egsl_atmp(eigenvectors[j],i,0) = gsl_matrix_get(evec,i,j);
158  }
159 
160 
161  gsl_vector_free(eval);
162  gsl_matrix_free(evec);
163  gsl_matrix_free(A);
164 }
165 
166 
167 
168 
169 
170 
171 
172 
val egsl_inverse(val v1)
Definition: egsl_ops.c:115
void egsl_symm_eig(val v, double *eigenvalues, val *eigenvectors)
Definition: egsl_ops.c:133
#define m3(v1, v2, v3)
Definition: egsl_macros.h:14
val egsl_copy_val(val v1)
Definition: egsl_ops.c:69
Definition: egsl.h:12
void egsl_add_to_col(val v1, size_t j, val v2)
Definition: egsl_ops.c:54
val egsl_alloc(size_t rows, size_t columns)
Definition: egsl.c:159
#define m(v1, v2)
Definition: egsl_macros.h:13
val egsl_transpose(val v1)
Definition: egsl_ops.c:107
val egsl_mult(val v1, val v2)
Definition: egsl_ops.c:98
val egsl_sum3(val v1, val v2, val v3)
Definition: egsl_ops.c:94
val egsl_compose_col(val v1, val v2)
Definition: egsl_ops.c:14
void egsl_expect_size(val v, size_t rows, size_t cols)
Definition: egsl.c:238
gsl_matrix * egsl_gslm(val v)
Definition: egsl.c:83
double * egsl_atmp(val v, size_t i, size_t j)
Definition: egsl.c:276
val egsl_sum(val v1, val v2)
Definition: egsl_ops.c:84
val egsl_sub(val v1, val v2)
Definition: egsl_ops.c:10
val egsl_scale(double s, val v1)
Definition: egsl_ops.c:77
val egsl_compose_row(val v1, val v2)
Definition: egsl_ops.c:31
void egsl_add_to(val v1, val v2)
Definition: egsl_ops.c:48


csm
Author(s): Andrea Censi
autogenerated on Tue May 11 2021 02:18:23