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 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167
|
/* linalg/test_tri.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_blas.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_rng.h>
static int
test_symmtd_decomp_eps(const gsl_matrix * m, const double eps, const char * desc)
{
int s = 0;
const size_t N = m->size1;
size_t i, j;
gsl_matrix * Q = gsl_matrix_alloc(N, N);
gsl_matrix * T = gsl_matrix_calloc(N, N);
gsl_matrix * A = gsl_matrix_alloc(N, N);
gsl_matrix * B = gsl_matrix_alloc(N, N);
gsl_vector * tau = gsl_vector_alloc(N - 1);
gsl_vector_view diag = gsl_matrix_diagonal(T);
gsl_vector_view subdiag = gsl_matrix_subdiagonal(T, 1);
gsl_vector_view superdiag = gsl_matrix_superdiagonal(T, 1);
gsl_matrix_memcpy(A, m);
s += gsl_linalg_symmtd_decomp(A, tau);
s += gsl_linalg_symmtd_unpack(A, tau, Q, &diag.vector, &subdiag.vector);
gsl_vector_memcpy(&superdiag.vector, &subdiag.vector);
gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, Q, T, 0.0, A); /* A := Q T */
gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, A, Q, 0.0, B); /* B := Q T Q^T */
for (i = 0; i < N; i++)
{
for (j = 0; j <= i; j++)
{
double bij = gsl_matrix_get(B, i, j);
double mij = gsl_matrix_get(m, i, j);
gsl_test_rel(bij, mij, eps, "%s (%3lu)[%lu,%lu]: %22.18g %22.18g\n",
desc, N, i, j, bij, mij);
}
}
gsl_matrix_free(T);
gsl_matrix_free(A);
gsl_matrix_free(Q);
gsl_matrix_free(B);
gsl_vector_free(tau);
return s;
}
static int
test_symmtd_decomp(gsl_rng * r)
{
int s = 0;
size_t N;
for (N = 2; N <= 50; ++N)
{
gsl_matrix * A = gsl_matrix_alloc(N, N);
create_symm_matrix(A, r);
s += test_symmtd_decomp_eps(A, 1.0e5 * N * GSL_DBL_EPSILON, "symmtd_decomp random");
gsl_matrix_free(A);
}
return s;
}
static int
test_hermtd_decomp_eps(const gsl_matrix_complex * m, const double eps, const char * desc)
{
int s = 0;
const size_t N = m->size1;
size_t i, j;
gsl_matrix_complex * Q = gsl_matrix_complex_alloc(N, N);
gsl_matrix_complex * T = gsl_matrix_complex_calloc(N, N);
gsl_matrix_complex * A = gsl_matrix_complex_alloc(N, N);
gsl_matrix_complex * B = gsl_matrix_complex_alloc(N, N);
gsl_vector_complex * tau = gsl_vector_complex_alloc(N - 1);
gsl_vector_view diag, subdiag, superdiag;
gsl_vector_complex_view v;
v = gsl_matrix_complex_diagonal(T);
diag = gsl_vector_complex_real(&v.vector);
v = gsl_matrix_complex_subdiagonal(T, 1);
subdiag = gsl_vector_complex_real(&v.vector);
v = gsl_matrix_complex_superdiagonal(T, 1);
superdiag = gsl_vector_complex_real(&v.vector);
gsl_matrix_complex_memcpy(A, m);
s += gsl_linalg_hermtd_decomp(A, tau);
s += gsl_linalg_hermtd_unpack(A, tau, Q, &diag.vector, &subdiag.vector);
gsl_vector_memcpy(&superdiag.vector, &subdiag.vector);
gsl_blas_zgemm(CblasNoTrans, CblasNoTrans, GSL_COMPLEX_ONE, Q, T, GSL_COMPLEX_ZERO, A); /* A := Q T */
gsl_blas_zgemm(CblasNoTrans, CblasConjTrans, GSL_COMPLEX_ONE, A, Q, GSL_COMPLEX_ZERO, B); /* B := Q T Q^T */
for (i = 0; i < N; i++)
{
for (j = 0; j <= i; j++)
{
gsl_complex bij = gsl_matrix_complex_get(B, i, j);
gsl_complex mij = gsl_matrix_complex_get(m, i, j);
gsl_test_rel(GSL_REAL(bij), GSL_REAL(mij), eps, "%s real (%3lu)[%lu,%lu]: %22.18g %22.18g\n",
desc, N, i, j, GSL_REAL(bij), GSL_REAL(mij));
gsl_test_rel(GSL_IMAG(bij), GSL_IMAG(mij), eps, "%s imag (%3lu)[%lu,%lu]: %22.18g %22.18g\n",
desc, N, i, j, GSL_IMAG(bij), GSL_IMAG(mij));
}
}
gsl_matrix_complex_free(T);
gsl_matrix_complex_free(A);
gsl_matrix_complex_free(Q);
gsl_matrix_complex_free(B);
gsl_vector_complex_free(tau);
return s;
}
static int
test_hermtd_decomp(gsl_rng * r)
{
int s = 0;
size_t N;
for (N = 2; N <= 50; ++N)
{
gsl_matrix_complex * A = gsl_matrix_complex_alloc(N, N);
create_herm_matrix(A, r);
s += test_hermtd_decomp_eps(A, 1.0e5 * N * GSL_DBL_EPSILON, "hermtd_decomp random");
gsl_matrix_complex_free(A);
}
return s;
}
|