1 /* linalg/test_common.c
2  *
3  * Copyright (C) 2017, 2018, 2019, 2020 Patrick Alken
4  *
5  * This program is free software; you can redistribute it and/or modify
6  * it under the terms of the GNU General Public License as published by
7  * the Free Software Foundation; either version 3 of the License, or (at
8  * your option) any later version.
9  *
10  * This program is distributed in the hope that it will be useful, but
11  * WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
13  * General Public License for more details.
14  *
15  * You should have received a copy of the GNU General Public License
16  * along with this program; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
18  */
19 
20 #include <config.h>
21 #include <stdlib.h>
22 #include <gsl/gsl_math.h>
23 #include <gsl/gsl_blas.h>
24 #include <gsl/gsl_matrix.h>
25 #include <gsl/gsl_vector.h>
26 #include <gsl/gsl_rng.h>
27 
28 static int create_random_vector(gsl_vector * v, gsl_rng * r);
29 static int create_random_matrix(gsl_matrix * m, gsl_rng * r);
30 static int create_posdef_matrix(gsl_matrix * m, gsl_rng * r);
31 static int create_hilbert_matrix2(gsl_matrix * m);
32 
33 static int
create_random_vector(gsl_vector * v,gsl_rng * r)34 create_random_vector(gsl_vector * v, gsl_rng * r)
35 {
36   const size_t N = v->size;
37   size_t i;
38 
39   for (i = 0; i < N; ++i)
40     {
41       double vi = gsl_rng_uniform(r);
42       gsl_vector_set(v, i, vi);
43     }
44 
45   return GSL_SUCCESS;
46 }
47 
48 static int
create_random_complex_vector(gsl_vector_complex * v,gsl_rng * r)49 create_random_complex_vector(gsl_vector_complex * v, gsl_rng * r)
50 {
51   const size_t N = v->size;
52   size_t i;
53 
54   for (i = 0; i < N; ++i)
55     {
56       gsl_complex vi;
57       GSL_REAL(vi) = gsl_rng_uniform(r);
58       GSL_IMAG(vi) = gsl_rng_uniform(r);
59       gsl_vector_complex_set(v, i, vi);
60     }
61 
62   return GSL_SUCCESS;
63 }
64 
65 static int
create_random_matrix(gsl_matrix * m,gsl_rng * r)66 create_random_matrix(gsl_matrix * m, gsl_rng * r)
67 {
68   const size_t M = m->size1;
69   const size_t N = m->size2;
70   size_t i, j;
71 
72   for (i = 0; i < M; ++i)
73     {
74       for (j = 0; j < N; ++j)
75         {
76           double mij = gsl_rng_uniform(r);
77           gsl_matrix_set(m, i, j, mij);
78         }
79     }
80 
81   return GSL_SUCCESS;
82 }
83 
84 static int
create_random_complex_matrix(gsl_matrix_complex * m,gsl_rng * r)85 create_random_complex_matrix(gsl_matrix_complex * m, gsl_rng * r)
86 {
87   const size_t M = m->size1;
88   const size_t N = m->size2;
89   size_t i, j;
90 
91   for (i = 0; i < M; ++i)
92     {
93       for (j = 0; j < N; ++j)
94         {
95           gsl_complex mij;
96 
97           GSL_REAL(mij) = gsl_rng_uniform(r);
98           GSL_IMAG(mij) = gsl_rng_uniform(r);
99 
100           gsl_matrix_complex_set(m, i, j, mij);
101         }
102     }
103 
104   return GSL_SUCCESS;
105 }
106 
107 static int
create_symm_matrix(gsl_matrix * m,gsl_rng * r)108 create_symm_matrix(gsl_matrix * m, gsl_rng * r)
109 {
110   const size_t N = m->size1;
111   size_t i, j;
112 
113   for (i = 0; i < N; ++i)
114     {
115       for (j = 0; j <= i; ++j)
116         {
117           double mij = gsl_rng_uniform(r);
118           gsl_matrix_set(m, i, j, mij);
119         }
120     }
121 
122   /* copy lower triangle to upper */
123   gsl_matrix_transpose_tricpy(CblasLower, CblasUnit, m, m);
124 
125   return GSL_SUCCESS;
126 }
127 
128 static int
create_herm_matrix(gsl_matrix_complex * m,gsl_rng * r)129 create_herm_matrix(gsl_matrix_complex * m, gsl_rng * r)
130 {
131   const size_t N = m->size1;
132   size_t i, j;
133 
134   for (i = 0; i < N; ++i)
135     {
136       for (j = 0; j <= i; ++j)
137         {
138           double re = gsl_rng_uniform(r);
139           double im = (i != j) ? gsl_rng_uniform(r) : 0.0;
140           gsl_complex z = gsl_complex_rect(re, im);
141 
142           gsl_matrix_complex_set(m, i, j, z);
143 
144           if (i != j)
145             gsl_matrix_complex_set(m, j, i, gsl_complex_conjugate(z));
146         }
147     }
148 
149   return GSL_SUCCESS;
150 }
151 
152 /* create symmetric banded matrix with p sub/super-diagonals */
153 static int
create_symm_band_matrix(const size_t p,gsl_matrix * m,gsl_rng * r)154 create_symm_band_matrix(const size_t p, gsl_matrix * m, gsl_rng * r)
155 {
156   size_t i;
157 
158   gsl_matrix_set_zero(m);
159 
160   for (i = 0; i < p + 1; ++i)
161     {
162       gsl_vector_view subdiag = gsl_matrix_subdiagonal(m, i);
163       create_random_vector(&subdiag.vector, r);
164 
165       if (i > 0)
166         {
167           gsl_vector_view superdiag = gsl_matrix_superdiagonal(m, i);
168           gsl_vector_memcpy(&superdiag.vector, &subdiag.vector);
169         }
170     }
171 
172   return GSL_SUCCESS;
173 }
174 
175 /* create (p,q) banded matrix */
176 static int
create_band_matrix(const size_t p,const size_t q,gsl_matrix * m,gsl_rng * r)177 create_band_matrix(const size_t p, const size_t q, gsl_matrix * m, gsl_rng * r)
178 {
179   size_t i;
180 
181   gsl_matrix_set_zero(m);
182 
183   for (i = 0; i <= p; ++i)
184     {
185       gsl_vector_view v = gsl_matrix_subdiagonal(m, i);
186       create_random_vector(&v.vector, r);
187     }
188 
189   for (i = 1; i <= q; ++i)
190     {
191       gsl_vector_view v = gsl_matrix_superdiagonal(m, i);
192       create_random_vector(&v.vector, r);
193     }
194 
195   return GSL_SUCCESS;
196 }
197 
198 static int
create_posdef_matrix(gsl_matrix * m,gsl_rng * r)199 create_posdef_matrix(gsl_matrix * m, gsl_rng * r)
200 {
201   const size_t N = m->size1;
202   const double alpha = 10.0 * N;
203   size_t i;
204 
205   /* The idea is to make a symmetric diagonally dominant
206    * matrix. Make a symmetric matrix and add alpha*I to
207    * its diagonal */
208 
209   create_symm_matrix(m, r);
210 
211   for (i = 0; i < N; ++i)
212     {
213       double mii = gsl_matrix_get(m, i, i);
214       gsl_matrix_set(m, i, i, mii + alpha);
215     }
216 
217   return GSL_SUCCESS;
218 }
219 
220 static int
create_posdef_complex_matrix(gsl_matrix_complex * m,gsl_rng * r)221 create_posdef_complex_matrix(gsl_matrix_complex *m, gsl_rng *r)
222 {
223   const size_t N = m->size1;
224   const double alpha = 10.0 * N;
225   size_t i;
226 
227   create_herm_matrix(m, r);
228 
229   for (i = 0; i < N; ++i)
230     {
231       gsl_complex * mii = gsl_matrix_complex_ptr(m, i, i);
232       GSL_REAL(*mii) += alpha;
233     }
234 
235   return GSL_SUCCESS;
236 }
237 
238 static int
create_hilbert_matrix2(gsl_matrix * m)239 create_hilbert_matrix2(gsl_matrix * m)
240 {
241   const size_t N = m->size1;
242   size_t i, j;
243 
244   for (i = 0; i < N; i++)
245     {
246       for (j = 0; j < N; j++)
247         {
248           gsl_matrix_set(m, i, j, 1.0/(i+j+1.0));
249         }
250     }
251 
252   return GSL_SUCCESS;
253 }
254 
255 static int
create_posdef_band_matrix(const size_t p,gsl_matrix * m,gsl_rng * r)256 create_posdef_band_matrix(const size_t p, gsl_matrix * m, gsl_rng * r)
257 {
258   const size_t N = m->size1;
259   const double alpha = 10.0 * N;
260   size_t i;
261 
262   /* The idea is to make a symmetric diagonally dominant
263    * matrix. Make a symmetric matrix and add alpha*I to
264    * its diagonal */
265 
266   create_symm_band_matrix(p, m, r);
267 
268   for (i = 0; i < N; ++i)
269     {
270       double *mii = gsl_matrix_ptr(m, i, i);
271       *mii += alpha;
272     }
273 
274   return GSL_SUCCESS;
275 }
276 
277 /* transform dense symmetric banded matrix to compact form, with bandwidth p */
278 static int
symm2band_matrix(const size_t p,const gsl_matrix * m,gsl_matrix * bm)279 symm2band_matrix(const size_t p, const gsl_matrix * m, gsl_matrix * bm)
280 {
281   const size_t N = m->size1;
282 
283   if (bm->size1 != N)
284     {
285       GSL_ERROR("banded matrix requires N rows", GSL_EBADLEN);
286     }
287   else if (bm->size2 != p + 1)
288     {
289       GSL_ERROR("banded matrix requires p + 1 columns", GSL_EBADLEN);
290     }
291   else
292     {
293       size_t i;
294 
295       gsl_matrix_set_zero(bm);
296 
297       for (i = 0; i < p + 1; ++i)
298         {
299           gsl_vector_const_view diag = gsl_matrix_const_subdiagonal(m, i);
300           gsl_vector_view v = gsl_matrix_subcolumn(bm, i, 0, N - i);
301 
302           gsl_vector_memcpy(&v.vector, &diag.vector);
303         }
304 
305       return GSL_SUCCESS;
306     }
307 }
308 
309 /* transform general dense (p,q) banded matrix to compact banded form */
310 static int
gen2band_matrix(const size_t p,const size_t q,const gsl_matrix * A,gsl_matrix * AB)311 gen2band_matrix(const size_t p, const size_t q, const gsl_matrix * A, gsl_matrix * AB)
312 {
313   const size_t N = A->size2;
314 
315   if (AB->size1 != N)
316     {
317       GSL_ERROR("banded matrix requires N rows", GSL_EBADLEN);
318     }
319   else if (AB->size2 != 2*p + q + 1)
320     {
321       GSL_ERROR("banded matrix requires 2*p + q + 1 columns", GSL_EBADLEN);
322     }
323   else
324     {
325       size_t i;
326 
327       gsl_matrix_set_zero(AB);
328 
329       /* copy diagonal and subdiagonals */
330       for (i = 0; i <= p; ++i)
331         {
332           gsl_vector_const_view v = gsl_matrix_const_subdiagonal(A, i);
333           gsl_vector_view w = gsl_matrix_subcolumn(AB, p + q + i, 0, v.vector.size);
334           gsl_vector_memcpy(&w.vector, &v.vector);
335         }
336 
337       /* copy superdiagonals */
338       for (i = 1; i <= q; ++i)
339         {
340           gsl_vector_const_view v = gsl_matrix_const_superdiagonal(A, i);
341           gsl_vector_view w = gsl_matrix_subcolumn(AB, p + q - i, i, v.vector.size);
342           gsl_vector_memcpy(&w.vector, &v.vector);
343         }
344 
345       return GSL_SUCCESS;
346     }
347 }
348