Blame eigen/genhermv.c

Packit 67cb25
/* eigen/genhermv.c
Packit 67cb25
 * 
Packit 67cb25
 * Copyright (C) 2007 Patrick Alken
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
#include <stdlib.h>
Packit 67cb25
Packit 67cb25
#include <config.h>
Packit 67cb25
#include <gsl/gsl_eigen.h>
Packit 67cb25
#include <gsl/gsl_linalg.h>
Packit 67cb25
#include <gsl/gsl_math.h>
Packit 67cb25
#include <gsl/gsl_blas.h>
Packit 67cb25
#include <gsl/gsl_vector.h>
Packit 67cb25
#include <gsl/gsl_matrix.h>
Packit 67cb25
#include <gsl/gsl_complex.h>
Packit 67cb25
#include <gsl/gsl_complex_math.h>
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
 * This module computes the eigenvalues and eigenvectors of a complex
Packit 67cb25
 * generalized hermitian-definite eigensystem A x = \lambda B x, where
Packit 67cb25
 * A and B are hermitian, and B is positive-definite.
Packit 67cb25
 */
Packit 67cb25
Packit 67cb25
static void genhermv_normalize_eigenvectors(gsl_matrix_complex *evec);
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
gsl_eigen_genhermv_alloc()
Packit 67cb25
Packit 67cb25
Allocate a workspace for solving the generalized hermitian-definite
Packit 67cb25
eigenvalue problem. The size of this workspace is O(5n).
Packit 67cb25
Packit 67cb25
Inputs: n - size of matrices
Packit 67cb25
Packit 67cb25
Return: pointer to workspace
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
gsl_eigen_genhermv_workspace *
Packit 67cb25
gsl_eigen_genhermv_alloc(const size_t n)
Packit 67cb25
{
Packit 67cb25
  gsl_eigen_genhermv_workspace *w;
Packit 67cb25
Packit 67cb25
  if (n == 0)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL ("matrix dimension must be positive integer",
Packit 67cb25
                      GSL_EINVAL);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  w = (gsl_eigen_genhermv_workspace *) calloc (1, sizeof (gsl_eigen_genhermv_workspace));
Packit 67cb25
Packit 67cb25
  if (w == 0)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR_NULL ("failed to allocate space for workspace", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  w->size = n;
Packit 67cb25
Packit 67cb25
  w->hermv_workspace_p = gsl_eigen_hermv_alloc(n);
Packit 67cb25
  if (!w->hermv_workspace_p)
Packit 67cb25
    {
Packit 67cb25
      gsl_eigen_genhermv_free(w);
Packit 67cb25
      GSL_ERROR_NULL("failed to allocate space for hermv workspace", GSL_ENOMEM);
Packit 67cb25
    }
Packit 67cb25
Packit 67cb25
  return (w);
Packit 67cb25
} /* gsl_eigen_genhermv_alloc() */
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
gsl_eigen_genhermv_free()
Packit 67cb25
  Free workspace w
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
void
Packit 67cb25
gsl_eigen_genhermv_free (gsl_eigen_genhermv_workspace * w)
Packit 67cb25
{
Packit 67cb25
  RETURN_IF_NULL (w);
Packit 67cb25
Packit 67cb25
  if (w->hermv_workspace_p)
Packit 67cb25
    gsl_eigen_hermv_free(w->hermv_workspace_p);
Packit 67cb25
Packit 67cb25
  free(w);
Packit 67cb25
} /* gsl_eigen_genhermv_free() */
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
gsl_eigen_genhermv()
Packit 67cb25
Packit 67cb25
Solve the generalized hermitian-definite eigenvalue problem
Packit 67cb25
Packit 67cb25
A x = \lambda B x
Packit 67cb25
Packit 67cb25
for the eigenvalues \lambda and eigenvectors x.
Packit 67cb25
Packit 67cb25
Inputs: A    - complex hermitian matrix
Packit 67cb25
        B    - complex hermitian and positive definite matrix
Packit 67cb25
        eval - where to store eigenvalues
Packit 67cb25
        evec - where to store eigenvectors
Packit 67cb25
        w    - workspace
Packit 67cb25
Packit 67cb25
Return: success or error
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
int
Packit 67cb25
gsl_eigen_genhermv (gsl_matrix_complex * A, gsl_matrix_complex * B,
Packit 67cb25
                    gsl_vector * eval, gsl_matrix_complex * evec,
Packit 67cb25
                    gsl_eigen_genhermv_workspace * w)
Packit 67cb25
{
Packit 67cb25
  const size_t N = A->size1;
Packit 67cb25
Packit 67cb25
  /* check matrix and vector sizes */
Packit 67cb25
Packit 67cb25
  if (N != A->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix must be square to compute eigenvalues", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if ((N != B->size1) || (N != B->size2))
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("B matrix dimensions must match A", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (eval->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("eigenvalue vector must match matrix size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (evec->size1 != evec->size2)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("eigenvector matrix must be square", GSL_ENOTSQR);
Packit 67cb25
    }
Packit 67cb25
  else if (evec->size1 != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("eigenvector matrix has wrong size", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else if (w->size != N)
Packit 67cb25
    {
Packit 67cb25
      GSL_ERROR ("matrix size does not match workspace", GSL_EBADLEN);
Packit 67cb25
    }
Packit 67cb25
  else
Packit 67cb25
    {
Packit 67cb25
      int s;
Packit 67cb25
Packit 67cb25
      /* compute Cholesky factorization of B */
Packit 67cb25
      s = gsl_linalg_complex_cholesky_decomp(B);
Packit 67cb25
      if (s != GSL_SUCCESS)
Packit 67cb25
        return s; /* B is not positive definite */
Packit 67cb25
Packit 67cb25
      /* transform to standard hermitian eigenvalue problem */
Packit 67cb25
      gsl_eigen_genherm_standardize(A, B);
Packit 67cb25
Packit 67cb25
      /* compute eigenvalues and eigenvectors */
Packit 67cb25
      s = gsl_eigen_hermv(A, eval, evec, w->hermv_workspace_p);
Packit 67cb25
      if (s != GSL_SUCCESS)
Packit 67cb25
        return s;
Packit 67cb25
Packit 67cb25
      /* backtransform eigenvectors: evec -> L^{-H} evec */
Packit 67cb25
      gsl_blas_ztrsm(CblasLeft,
Packit 67cb25
                     CblasLower,
Packit 67cb25
                     CblasConjTrans,
Packit 67cb25
                     CblasNonUnit,
Packit 67cb25
                     GSL_COMPLEX_ONE,
Packit 67cb25
                     B,
Packit 67cb25
                     evec);
Packit 67cb25
Packit 67cb25
      /* the blas call destroyed the normalization - renormalize */
Packit 67cb25
      genhermv_normalize_eigenvectors(evec);
Packit 67cb25
Packit 67cb25
      return GSL_SUCCESS;
Packit 67cb25
    }
Packit 67cb25
} /* gsl_eigen_genhermv() */
Packit 67cb25
Packit 67cb25
/********************************************
Packit 67cb25
 *           INTERNAL ROUTINES              *
Packit 67cb25
 ********************************************/
Packit 67cb25
Packit 67cb25
/*
Packit 67cb25
genhermv_normalize_eigenvectors()
Packit 67cb25
  Normalize eigenvectors so that their Euclidean norm is 1
Packit 67cb25
Packit 67cb25
Inputs: evec - eigenvectors
Packit 67cb25
*/
Packit 67cb25
Packit 67cb25
static void
Packit 67cb25
genhermv_normalize_eigenvectors(gsl_matrix_complex *evec)
Packit 67cb25
{
Packit 67cb25
  const size_t N = evec->size1;
Packit 67cb25
  size_t i;     /* looping */
Packit 67cb25
Packit 67cb25
  for (i = 0; i < N; ++i)
Packit 67cb25
    {
Packit 67cb25
      gsl_vector_complex_view vi = gsl_matrix_complex_column(evec, i);
Packit 67cb25
      double scale = 1.0 / gsl_blas_dznrm2(&vi.vector);
Packit 67cb25
Packit 67cb25
      gsl_blas_zdscal(scale, &vi.vector);
Packit 67cb25
    }
Packit 67cb25
} /* genhermv_normalize_eigenvectors() */