File: test_ql.c

package info (click to toggle)
gsl 2.7.1%2Bdfsg-5%2Bdeb12u1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 26,668 kB
  • sloc: ansic: 259,467; sh: 4,535; makefile: 902; python: 69
file content (103 lines) | stat: -rw-r--r-- 3,175 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
/* linalg/test_ql.c
 *
 * Copyright (C) 2019 Patrick Alken
 * 
 * This program is free software; you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation; either version 3 of the License, or (at
 * your option) any later version.
 * 
 * This program is distributed in the hope that it will be useful, but
 * WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
 * General Public License for more details.
 * 
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, write to the Free Software
 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
 */

#include <config.h>
#include <stdlib.h>
#include <gsl/gsl_test.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_ieee_utils.h>
#include <gsl/gsl_permute_vector.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_rng.h>
#include <gsl/gsl_permutation.h>

static int
test_QL_decomp_eps(const gsl_matrix * m, const double eps, const char * desc)
{
  int s = 0;
  const size_t M = m->size1;
  const size_t N = m->size2;
  size_t i, j;

  gsl_matrix * QL = gsl_matrix_alloc(M, N);
  gsl_vector * tau = gsl_vector_alloc(N);
  gsl_matrix * A  = gsl_matrix_alloc(M, N);
  gsl_matrix * L  = gsl_matrix_alloc(M, N);
  gsl_matrix * Q  = gsl_matrix_alloc(M, M);

  gsl_matrix_memcpy(QL, m);

  s += gsl_linalg_QL_decomp(QL, tau);
  s += gsl_linalg_QL_unpack(QL, tau, Q, L);

  /* compute A = Q L */
  gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0, Q, L, 0.0, A);

  for (i = 0; i < M; i++)
    {
      for (j = 0; j < N; j++)
        {
          double aij = gsl_matrix_get(A, i, j);
          double mij = gsl_matrix_get(m, i, j);

          gsl_test_rel(aij, mij, eps, "%s (%3lu,%3lu)[%lu,%lu]: %22.18g   %22.18g\n",
                       desc, M, N, i,j, aij, mij);
        }
    }

  gsl_matrix_free(QL);
  gsl_vector_free(tau);
  gsl_matrix_free(A);
  gsl_matrix_free(Q);
  gsl_matrix_free(L);

  return s;
}

static int
test_QL_decomp(gsl_rng * r)
{
  int s = 0;
  size_t M, N;

  for (M = 1; M <= 30; ++M)
    {
      for (N = 1; N <= M; ++N)
        {
          gsl_matrix * A = gsl_matrix_alloc(M, N);

          create_random_matrix(A, r);
          s += test_QL_decomp_eps(A, 1.0e5 * GSL_MAX(M, N) * GSL_DBL_EPSILON, "QL_decomp random");

          gsl_matrix_free(A);
        }
    }

  s += test_QL_decomp_eps(m53,   1.0e2 * GSL_DBL_EPSILON, "QL_decomp m(5,3)");
  s += test_QL_decomp_eps(hilb2, 1.0e2 * GSL_DBL_EPSILON, "QL_decomp hilbert(2)");
  s += test_QL_decomp_eps(hilb3, 1.0e2 * GSL_DBL_EPSILON, "QL_decomp hilbert(3)");
  s += test_QL_decomp_eps(hilb4, 1.0e2 * GSL_DBL_EPSILON, "QL_decomp hilbert(4)");
  s += test_QL_decomp_eps(hilb12, 1.0e2 * GSL_DBL_EPSILON, "QL_decomp hilbert(12)");
  s += test_QL_decomp_eps(vander2, 1.0e1 * GSL_DBL_EPSILON, "QL_decomp vander(2)");
  s += test_QL_decomp_eps(vander3, 1.0e1 * GSL_DBL_EPSILON, "QL_decomp vander(3)");
  s += test_QL_decomp_eps(vander4, 1.0e2 * GSL_DBL_EPSILON, "QL_decomp vander(4)");

  return s;
}