This is the mail archive of the gsl-discuss@sources.redhat.com mailing list for the GSL project.


Index Nav: [Date Index] [Subject Index] [Author Index] [Thread Index]
Message Nav: [Date Prev] [Date Next] [Thread Prev] [Thread Next]
Other format: [Raw text]

again gsl_linalg_complex_householder_transform


Hi,
 changing in gsl_linalg_complex_householder_transform
 GSL_IMAG(tau) = - GSL_IMAG(alpha) / beta_r ;
 into
 GSL_IMAG(tau) = GSL_IMAG(alpha) / beta_r ;
 the following example works, so forget about my previous modification
 and comments on this function.
 Forget also about my remarks on gsl_linalg_hermtd_decomp;
 it should give a real symmetric tridiagonal matrix, but I
 can't make yet a working example :(
 Bye
  Mario

#include <stdio.h>
#include <gsl/gsl_math.h> 
#include <gsl/gsl_vector.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_blas.h>
#include <gsl/gsl_complex_math.h>
#include <gsl/gsl_linalg.h>

int main()
{
 gsl_complex zero, one, c1, c2;
 gsl_vector_complex * v, * vv, * res;
 gsl_complex tau, cv0, cv1;
zero = gsl_complex_rect(0.0, 0.0);
 one = gsl_complex_rect(1.0, 0.0);
 c1 = gsl_complex_rect(3, -2);
 c2 = gsl_complex_rect(4, -1);

 v = gsl_vector_complex_alloc(2);
 vv = gsl_vector_complex_alloc(2);
 res = gsl_vector_complex_alloc(2);
 gsl_vector_complex_set(v,0,c1);
 gsl_vector_complex_set(v,1,c2);
 gsl_vector_complex_memcpy(vv,v);
 gsl_vector_complex_memcpy(res,v);

 tau = gsl_linalg_complex_householder_transform(vv);
 puts("Householder vector");
 gsl_vector_complex_set(vv,0, one);
 gsl_vector_complex_fprintf(stdout, vv, "%g");
 puts("res");
 gsl_linalg_complex_householder_hv(tau, vv, res);
 gsl_vector_complex_fprintf(stdout, res, "%g");
}
==================
Output:
Householder vector
1 0
0.473337 -0.00629059
res
-5.47723 0
1.77636e-15 -2.22045e-16



Index Nav: [Date Index] [Subject Index] [Author Index] [Thread Index]
Message Nav: [Date Prev] [Date Next] [Thread Prev] [Thread Next]