|
Packit |
67cb25 |
/* linalg/tridiag.c
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* Copyright (C) 1996, 1997, 1998, 1999, 2000, 2002, 2004, 2007 Gerard Jungman, Brian Gough, David Necas
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* This program is free software; you can redistribute it and/or modify
|
|
Packit |
67cb25 |
* it under the terms of the GNU General Public License as published by
|
|
Packit |
67cb25 |
* the Free Software Foundation; either version 3 of the License, or (at
|
|
Packit |
67cb25 |
* your option) any later version.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* This program is distributed in the hope that it will be useful, but
|
|
Packit |
67cb25 |
* WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
Packit |
67cb25 |
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
|
|
Packit |
67cb25 |
* General Public License for more details.
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* You should have received a copy of the GNU General Public License
|
|
Packit |
67cb25 |
* along with this program; if not, write to the Free Software
|
|
Packit |
67cb25 |
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Author: G. Jungman */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
#include <config.h>
|
|
Packit |
67cb25 |
#include <stdlib.h>
|
|
Packit |
67cb25 |
#include <math.h>
|
|
Packit |
67cb25 |
#include <gsl/gsl_errno.h>
|
|
Packit |
67cb25 |
#include "tridiag.h"
|
|
Packit |
67cb25 |
#include <gsl/gsl_linalg.h>
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* for description of method see [Engeln-Mullges + Uhlig, p. 92]
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* diag[0] offdiag[0] 0 .....
|
|
Packit |
67cb25 |
* offdiag[0] diag[1] offdiag[1] .....
|
|
Packit |
67cb25 |
* 0 offdiag[1] diag[2]
|
|
Packit |
67cb25 |
* 0 0 offdiag[2] .....
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
static
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
solve_tridiag(
|
|
Packit |
67cb25 |
const double diag[], size_t d_stride,
|
|
Packit |
67cb25 |
const double offdiag[], size_t o_stride,
|
|
Packit |
67cb25 |
const double b[], size_t b_stride,
|
|
Packit |
67cb25 |
double x[], size_t x_stride,
|
|
Packit |
67cb25 |
size_t N)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = GSL_SUCCESS;
|
|
Packit |
67cb25 |
double *gamma = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *alpha = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *c = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *z = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (gamma == 0 || alpha == 0 || c == 0 || z == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i, j;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Cholesky decomposition
|
|
Packit |
67cb25 |
A = L.D.L^t
|
|
Packit |
67cb25 |
lower_diag(L) = gamma
|
|
Packit |
67cb25 |
diag(D) = alpha
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
alpha[0] = diag[0];
|
|
Packit |
67cb25 |
gamma[0] = offdiag[0] / alpha[0];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (alpha[0] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 1; i < N - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
alpha[i] = diag[d_stride * i] - offdiag[o_stride*(i - 1)] * gamma[i - 1];
|
|
Packit |
67cb25 |
gamma[i] = offdiag[o_stride * i] / alpha[i];
|
|
Packit |
67cb25 |
if (alpha[i] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (N > 1)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
alpha[N - 1] = diag[d_stride * (N - 1)] - offdiag[o_stride*(N - 2)] * gamma[N - 2];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* update RHS */
|
|
Packit |
67cb25 |
z[0] = b[0];
|
|
Packit |
67cb25 |
for (i = 1; i < N; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
z[i] = b[b_stride * i] - gamma[i - 1] * z[i - 1];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
for (i = 0; i < N; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
c[i] = z[i] / alpha[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* backsubstitution */
|
|
Packit |
67cb25 |
x[x_stride * (N - 1)] = c[N - 1];
|
|
Packit |
67cb25 |
if (N >= 2)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
for (i = N - 2, j = 0; j <= N - 2; j++, i--)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
x[x_stride * i] = c[i] - gamma[i] * x[x_stride * (i + 1)];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (z != 0)
|
|
Packit |
67cb25 |
free (z);
|
|
Packit |
67cb25 |
if (c != 0)
|
|
Packit |
67cb25 |
free (c);
|
|
Packit |
67cb25 |
if (alpha != 0)
|
|
Packit |
67cb25 |
free (alpha);
|
|
Packit |
67cb25 |
if (gamma != 0)
|
|
Packit |
67cb25 |
free (gamma);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (status == GSL_EZERODIV) {
|
|
Packit |
67cb25 |
GSL_ERROR ("matrix must be positive definite", status);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* plain gauss elimination, only not bothering with the zeroes
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* diag[0] abovediag[0] 0 .....
|
|
Packit |
67cb25 |
* belowdiag[0] diag[1] abovediag[1] .....
|
|
Packit |
67cb25 |
* 0 belowdiag[1] diag[2]
|
|
Packit |
67cb25 |
* 0 0 belowdiag[2] .....
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
static
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
solve_tridiag_nonsym(
|
|
Packit |
67cb25 |
const double diag[], size_t d_stride,
|
|
Packit |
67cb25 |
const double abovediag[], size_t a_stride,
|
|
Packit |
67cb25 |
const double belowdiag[], size_t b_stride,
|
|
Packit |
67cb25 |
const double rhs[], size_t r_stride,
|
|
Packit |
67cb25 |
double x[], size_t x_stride,
|
|
Packit |
67cb25 |
size_t N)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = GSL_SUCCESS;
|
|
Packit |
67cb25 |
double *alpha = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *z = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (alpha == 0 || z == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i, j;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Bidiagonalization (eliminating belowdiag)
|
|
Packit |
67cb25 |
& rhs update
|
|
Packit |
67cb25 |
diag' = alpha
|
|
Packit |
67cb25 |
rhs' = z
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
alpha[0] = diag[0];
|
|
Packit |
67cb25 |
z[0] = rhs[0];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (alpha[0] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 1; i < N; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const double t = belowdiag[b_stride*(i - 1)]/alpha[i-1];
|
|
Packit |
67cb25 |
alpha[i] = diag[d_stride*i] - t*abovediag[a_stride*(i - 1)];
|
|
Packit |
67cb25 |
z[i] = rhs[r_stride*i] - t*z[i-1];
|
|
Packit |
67cb25 |
if (alpha[i] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* backsubstitution */
|
|
Packit |
67cb25 |
x[x_stride * (N - 1)] = z[N - 1]/alpha[N - 1];
|
|
Packit |
67cb25 |
if (N >= 2)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
for (i = N - 2, j = 0; j <= N - 2; j++, i--)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
x[x_stride * i] = (z[i] - abovediag[a_stride*i] * x[x_stride * (i + 1)])/alpha[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (z != 0)
|
|
Packit |
67cb25 |
free (z);
|
|
Packit |
67cb25 |
if (alpha != 0)
|
|
Packit |
67cb25 |
free (alpha);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (status == GSL_EZERODIV) {
|
|
Packit |
67cb25 |
GSL_ERROR ("matrix must be positive definite", status);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* for description of method see [Engeln-Mullges + Uhlig, p. 96]
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* diag[0] offdiag[0] 0 ..... offdiag[N-1]
|
|
Packit |
67cb25 |
* offdiag[0] diag[1] offdiag[1] .....
|
|
Packit |
67cb25 |
* 0 offdiag[1] diag[2]
|
|
Packit |
67cb25 |
* 0 0 offdiag[2] .....
|
|
Packit |
67cb25 |
* ... ...
|
|
Packit |
67cb25 |
* offdiag[N-1] ...
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
static
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
solve_cyc_tridiag(
|
|
Packit |
67cb25 |
const double diag[], size_t d_stride,
|
|
Packit |
67cb25 |
const double offdiag[], size_t o_stride,
|
|
Packit |
67cb25 |
const double b[], size_t b_stride,
|
|
Packit |
67cb25 |
double x[], size_t x_stride,
|
|
Packit |
67cb25 |
size_t N)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = GSL_SUCCESS;
|
|
Packit |
67cb25 |
double * delta = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double * gamma = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double * alpha = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double * c = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double * z = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (delta == 0 || gamma == 0 || alpha == 0 || c == 0 || z == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i, j;
|
|
Packit |
67cb25 |
double sum = 0.0;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* factor */
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (N == 1)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
x[0] = b[0] / diag[0];
|
|
Packit |
67cb25 |
free(delta);
|
|
Packit |
67cb25 |
free(gamma);
|
|
Packit |
67cb25 |
free(alpha);
|
|
Packit |
67cb25 |
free(c);
|
|
Packit |
67cb25 |
free(z);
|
|
Packit |
67cb25 |
return GSL_SUCCESS;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
alpha[0] = diag[0];
|
|
Packit |
67cb25 |
gamma[0] = offdiag[0] / alpha[0];
|
|
Packit |
67cb25 |
delta[0] = offdiag[o_stride * (N-1)] / alpha[0];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (alpha[0] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 1; i < N - 2; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
alpha[i] = diag[d_stride * i] - offdiag[o_stride * (i-1)] * gamma[i - 1];
|
|
Packit |
67cb25 |
gamma[i] = offdiag[o_stride * i] / alpha[i];
|
|
Packit |
67cb25 |
delta[i] = -delta[i - 1] * offdiag[o_stride * (i-1)] / alpha[i];
|
|
Packit |
67cb25 |
if (alpha[i] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
for (i = 0; i < N - 2; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
sum += alpha[i] * delta[i] * delta[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
alpha[N - 2] = diag[d_stride * (N - 2)] - offdiag[o_stride * (N - 3)] * gamma[N - 3];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
gamma[N - 2] = (offdiag[o_stride * (N - 2)] - offdiag[o_stride * (N - 3)] * delta[N - 3]) / alpha[N - 2];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
alpha[N - 1] = diag[d_stride * (N - 1)] - sum - alpha[(N - 2)] * gamma[N - 2] * gamma[N - 2];
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* update */
|
|
Packit |
67cb25 |
z[0] = b[0];
|
|
Packit |
67cb25 |
for (i = 1; i < N - 1; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
z[i] = b[b_stride * i] - z[i - 1] * gamma[i - 1];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
sum = 0.0;
|
|
Packit |
67cb25 |
for (i = 0; i < N - 2; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
sum += delta[i] * z[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
z[N - 1] = b[b_stride * (N - 1)] - sum - gamma[N - 2] * z[N - 2];
|
|
Packit |
67cb25 |
for (i = 0; i < N; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
c[i] = z[i] / alpha[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* backsubstitution */
|
|
Packit |
67cb25 |
x[x_stride * (N - 1)] = c[N - 1];
|
|
Packit |
67cb25 |
x[x_stride * (N - 2)] = c[N - 2] - gamma[N - 2] * x[x_stride * (N - 1)];
|
|
Packit |
67cb25 |
if (N >= 3)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
for (i = N - 3, j = 0; j <= N - 3; j++, i--)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
x[x_stride * i] = c[i] - gamma[i] * x[x_stride * (i + 1)] - delta[i] * x[x_stride * (N - 1)];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (z != 0)
|
|
Packit |
67cb25 |
free (z);
|
|
Packit |
67cb25 |
if (c != 0)
|
|
Packit |
67cb25 |
free (c);
|
|
Packit |
67cb25 |
if (alpha != 0)
|
|
Packit |
67cb25 |
free (alpha);
|
|
Packit |
67cb25 |
if (gamma != 0)
|
|
Packit |
67cb25 |
free (gamma);
|
|
Packit |
67cb25 |
if (delta != 0)
|
|
Packit |
67cb25 |
free (delta);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (status == GSL_EZERODIV) {
|
|
Packit |
67cb25 |
GSL_ERROR ("matrix must be positive definite", status);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* solve following system w/o the corner elements and then use
|
|
Packit |
67cb25 |
* Sherman-Morrison formula to compensate for them
|
|
Packit |
67cb25 |
*
|
|
Packit |
67cb25 |
* diag[0] abovediag[0] 0 ..... belowdiag[N-1]
|
|
Packit |
67cb25 |
* belowdiag[0] diag[1] abovediag[1] .....
|
|
Packit |
67cb25 |
* 0 belowdiag[1] diag[2]
|
|
Packit |
67cb25 |
* 0 0 belowdiag[2] .....
|
|
Packit |
67cb25 |
* ... ...
|
|
Packit |
67cb25 |
* abovediag[N-1] ...
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
static
|
|
Packit |
67cb25 |
int solve_cyc_tridiag_nonsym(
|
|
Packit |
67cb25 |
const double diag[], size_t d_stride,
|
|
Packit |
67cb25 |
const double abovediag[], size_t a_stride,
|
|
Packit |
67cb25 |
const double belowdiag[], size_t b_stride,
|
|
Packit |
67cb25 |
const double rhs[], size_t r_stride,
|
|
Packit |
67cb25 |
double x[], size_t x_stride,
|
|
Packit |
67cb25 |
size_t N)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
int status = GSL_SUCCESS;
|
|
Packit |
67cb25 |
double *alpha = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *zb = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *zu = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
double *w = (double *) malloc (N * sizeof (double));
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (alpha == 0 || zb == 0 || zu == 0 || w == 0)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR("failed to allocate working space", GSL_ENOMEM);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
double beta;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Bidiagonalization (eliminating belowdiag)
|
|
Packit |
67cb25 |
& rhs update
|
|
Packit |
67cb25 |
diag' = alpha
|
|
Packit |
67cb25 |
rhs' = zb
|
|
Packit |
67cb25 |
rhs' for Aq=u is zu
|
|
Packit |
67cb25 |
*/
|
|
Packit |
67cb25 |
zb[0] = rhs[0];
|
|
Packit |
67cb25 |
if (diag[0] != 0) beta = -diag[0]; else beta = 1;
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const double q = 1 - abovediag[0]*belowdiag[0]/(diag[0]*diag[d_stride]);
|
|
Packit |
67cb25 |
if (fabs(q/beta) > 0.5 && fabs(q/beta) < 2) {
|
|
Packit |
67cb25 |
beta *= (fabs(q/beta) < 1) ? 0.5 : 2;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
zu[0] = beta;
|
|
Packit |
67cb25 |
alpha[0] = diag[0] - beta;
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (alpha[0] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
for (i = 1; i+1 < N; i++)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const double t = belowdiag[b_stride*(i - 1)]/alpha[i-1];
|
|
Packit |
67cb25 |
alpha[i] = diag[d_stride*i] - t*abovediag[a_stride*(i - 1)];
|
|
Packit |
67cb25 |
zb[i] = rhs[r_stride*i] - t*zb[i-1];
|
|
Packit |
67cb25 |
zu[i] = -t*zu[i-1];
|
|
Packit |
67cb25 |
/* FIXME!!! */
|
|
Packit |
67cb25 |
if (alpha[i] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const size_t i = N-1;
|
|
Packit |
67cb25 |
const double t = belowdiag[b_stride*(i - 1)]/alpha[i-1];
|
|
Packit |
67cb25 |
alpha[i] = diag[d_stride*i]
|
|
Packit |
67cb25 |
- abovediag[a_stride*i]*belowdiag[b_stride*i]/beta
|
|
Packit |
67cb25 |
- t*abovediag[a_stride*(i - 1)];
|
|
Packit |
67cb25 |
zb[i] = rhs[r_stride*i] - t*zb[i-1];
|
|
Packit |
67cb25 |
zu[i] = abovediag[a_stride*i] - t*zu[i-1];
|
|
Packit |
67cb25 |
/* FIXME!!! */
|
|
Packit |
67cb25 |
if (alpha[i] == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* backsubstitution */
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i, j;
|
|
Packit |
67cb25 |
w[N-1] = zu[N-1]/alpha[N-1];
|
|
Packit |
67cb25 |
x[x_stride*(N-1)] = zb[N-1]/alpha[N-1];
|
|
Packit |
67cb25 |
for (i = N - 2, j = 0; j <= N - 2; j++, i--)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
w[i] = (zu[i] - abovediag[a_stride*i] * w[i+1])/alpha[i];
|
|
Packit |
67cb25 |
x[i*x_stride] = (zb[i] - abovediag[a_stride*i] * x[x_stride*(i + 1)])/alpha[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
/* Sherman-Morrison */
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
const double vw = w[0] + belowdiag[b_stride*(N - 1)]/beta * w[N-1];
|
|
Packit |
67cb25 |
const double vx = x[0] + belowdiag[b_stride*(N - 1)]/beta * x[x_stride*(N - 1)];
|
|
Packit |
67cb25 |
/* FIXME!!! */
|
|
Packit |
67cb25 |
if (vw + 1 == 0) {
|
|
Packit |
67cb25 |
status = GSL_EZERODIV;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
size_t i;
|
|
Packit |
67cb25 |
for (i = 0; i < N; i++)
|
|
Packit |
67cb25 |
x[i*x_stride] -= vx/(1 + vw)*w[i];
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (zb != 0)
|
|
Packit |
67cb25 |
free (zb);
|
|
Packit |
67cb25 |
if (zu != 0)
|
|
Packit |
67cb25 |
free (zu);
|
|
Packit |
67cb25 |
if (w != 0)
|
|
Packit |
67cb25 |
free (w);
|
|
Packit |
67cb25 |
if (alpha != 0)
|
|
Packit |
67cb25 |
free (alpha);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
if (status == GSL_EZERODIV) {
|
|
Packit |
67cb25 |
GSL_ERROR ("matrix must be positive definite", status);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
return status;
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_linalg_solve_symm_tridiag(
|
|
Packit |
67cb25 |
const gsl_vector * diag,
|
|
Packit |
67cb25 |
const gsl_vector * offdiag,
|
|
Packit |
67cb25 |
const gsl_vector * rhs,
|
|
Packit |
67cb25 |
gsl_vector * solution)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
if(diag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (offdiag->size != rhs->size-1)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of offdiag must match rhs-1", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (solution->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return solve_tridiag(diag->data, diag->stride,
|
|
Packit |
67cb25 |
offdiag->data, offdiag->stride,
|
|
Packit |
67cb25 |
rhs->data, rhs->stride,
|
|
Packit |
67cb25 |
solution->data, solution->stride,
|
|
Packit |
67cb25 |
diag->size);
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_linalg_solve_tridiag(
|
|
Packit |
67cb25 |
const gsl_vector * diag,
|
|
Packit |
67cb25 |
const gsl_vector * abovediag,
|
|
Packit |
67cb25 |
const gsl_vector * belowdiag,
|
|
Packit |
67cb25 |
const gsl_vector * rhs,
|
|
Packit |
67cb25 |
gsl_vector * solution)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
if(diag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (abovediag->size != rhs->size-1)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of abovediag must match rhs-1", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (belowdiag->size != rhs->size-1)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of belowdiag must match rhs-1", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (solution->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return solve_tridiag_nonsym(diag->data, diag->stride,
|
|
Packit |
67cb25 |
abovediag->data, abovediag->stride,
|
|
Packit |
67cb25 |
belowdiag->data, belowdiag->stride,
|
|
Packit |
67cb25 |
rhs->data, rhs->stride,
|
|
Packit |
67cb25 |
solution->data, solution->stride,
|
|
Packit |
67cb25 |
diag->size);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_linalg_solve_symm_cyc_tridiag(
|
|
Packit |
67cb25 |
const gsl_vector * diag,
|
|
Packit |
67cb25 |
const gsl_vector * offdiag,
|
|
Packit |
67cb25 |
const gsl_vector * rhs,
|
|
Packit |
67cb25 |
gsl_vector * solution)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
if(diag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (offdiag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of offdiag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (solution->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (diag->size < 3)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of cyclic system must be 3 or more", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return solve_cyc_tridiag(diag->data, diag->stride,
|
|
Packit |
67cb25 |
offdiag->data, offdiag->stride,
|
|
Packit |
67cb25 |
rhs->data, rhs->stride,
|
|
Packit |
67cb25 |
solution->data, solution->stride,
|
|
Packit |
67cb25 |
diag->size);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
int
|
|
Packit |
67cb25 |
gsl_linalg_solve_cyc_tridiag(
|
|
Packit |
67cb25 |
const gsl_vector * diag,
|
|
Packit |
67cb25 |
const gsl_vector * abovediag,
|
|
Packit |
67cb25 |
const gsl_vector * belowdiag,
|
|
Packit |
67cb25 |
const gsl_vector * rhs,
|
|
Packit |
67cb25 |
gsl_vector * solution)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
if(diag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of diag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (abovediag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of abovediag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (belowdiag->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of belowdiag must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (solution->size != rhs->size)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of solution must match rhs", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else if (diag->size < 3)
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
GSL_ERROR ("size of cyclic system must be 3 or more", GSL_EBADLEN);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
else
|
|
Packit |
67cb25 |
{
|
|
Packit |
67cb25 |
return solve_cyc_tridiag_nonsym(diag->data, diag->stride,
|
|
Packit |
67cb25 |
abovediag->data, abovediag->stride,
|
|
Packit |
67cb25 |
belowdiag->data, belowdiag->stride,
|
|
Packit |
67cb25 |
rhs->data, rhs->stride,
|
|
Packit |
67cb25 |
solution->data, solution->stride,
|
|
Packit |
67cb25 |
diag->size);
|
|
Packit |
67cb25 |
}
|
|
Packit |
67cb25 |
}
|