$include_dir="/home/hyper-archives/boost-commit/include"; include("$include_dir/msg-header.inc") ?>
From: thomas.klimpel_at_[hidden]
Date: 2008-05-15 15:00:58
Author: klimpel
Date: 2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
New Revision: 45398
URL: http://svn.boost.org/trac/boost/changeset/45398
Log:
added bindings for syevd, heevd, syevx and heevx
Added:
   sandbox/boost/numeric/bindings/lapack/heevd.hpp   (contents, props changed)
   sandbox/boost/numeric/bindings/lapack/heevx.hpp   (contents, props changed)
   sandbox/boost/numeric/bindings/lapack/syevd.hpp   (contents, props changed)
   sandbox/boost/numeric/bindings/lapack/syevx.hpp   (contents, props changed)
Text files modified: 
   sandbox/boost/numeric/bindings/lapack/lapack.h       |    46 ++++++++++++++++++++++++++++++++++++++++
   sandbox/boost/numeric/bindings/lapack/lapack_names.h |    10 ++++++++                                
   2 files changed, 56 insertions(+), 0 deletions(-)
Added: sandbox/boost/numeric/bindings/lapack/heevd.hpp
==============================================================================
--- (empty file)
+++ sandbox/boost/numeric/bindings/lapack/heevd.hpp	2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
@@ -0,0 +1,277 @@
+/*
+ * 
+ * Copyright Toon Knapen, Karl Meerbergen & Kresimir Fresl 2003
+ * Copyright Thomas Klimpel 2008
+ *
+ * Distributed under the Boost Software License, Version 1.0.
+ * (See accompanying file LICENSE_1_0.txt or copy at
+ * http://www.boost.org/LICENSE_1_0.txt)
+ *
+ * KF acknowledges the support of the Faculty of Civil Engineering,
+ * University of Zagreb, Croatia.
+ *
+ */
+
+#ifndef BOOST_NUMERIC_BINDINGS_LAPACK_HEEVD_HPP
+#define BOOST_NUMERIC_BINDINGS_LAPACK_HEEVD_HPP
+
+#include <boost/numeric/bindings/traits/traits.hpp>
+#include <boost/numeric/bindings/traits/type.hpp>
+#include <boost/numeric/bindings/lapack/lapack.h>
+#include <boost/numeric/bindings/lapack/workspace.hpp>
+#include <boost/numeric/bindings/traits/detail/array.hpp>
+
+#ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
+#  include <boost/static_assert.hpp>
+#  include <boost/type_traits.hpp>
+#endif
+
+
+namespace boost { namespace numeric { namespace bindings {
+
+  namespace lapack {
+
+    ///////////////////////////////////////////////////////////////////
+    //
+    // Eigendecomposition of a complex Hermitian matrix A = Q * D * Q'
+    //
+    ///////////////////////////////////////////////////////////////////
+
+    /*
+     * heevd() computes all eigenvalues and, optionally, eigenvectors of a
+     * complex Hermitian matrix A.  If eigenvectors are desired, it uses a
+     * divide and conquer algorithm.
+     *
+     * heevd() computes the eigendecomposition of a N x N matrix
+     * A = Q * D * Q',  where Q is a N x N unitary matrix and
+     * D is a diagonal matrix. The diagonal element D(i,i) is an
+     * eigenvalue of A and Q(:,i) is a corresponding eigenvector.
+     * The eigenvalues are stored in ascending order.
+     *
+     * On return of heevd, A is overwritten with the eigenvectors from Q
+     * and w contains the eigenvalues from the main diagonal of D.
+     *
+     * int heevd (char jobz, char uplo, A& a, W& w, Work work) {
+     *    jobz :  'V' : compute eigenvectors
+     *            'N' : do not compute eigenvectors
+     *    uplo :  'U' : only the upper triangular part of A is used on input.
+     *            'L' : only the lower triangular part of A is used on input.
+     */
+
+    namespace detail {
+
+      inline void heevd (
+        char const jobz, char const uplo, int const n,
+        float* a, int const lda,
+        float* w, float* work, int const lwork,
+        int* iwork, int const liwork, int& info)
+      {
+        LAPACK_SSYEVD (
+          &jobz, &uplo, &n,
+          a, &lda,
+          w, work, &lwork,
+          iwork, &liwork, &info);
+      }
+
+      inline void heevd (
+        char const jobz, char const uplo, int const n,
+        double* a, int const lda,
+        double* w, double* work, int const lwork,
+        int* iwork, int const liwork, int& info)
+      {
+        LAPACK_DSYEVD (
+          &jobz, &uplo, &n,
+          a, &lda,
+          w, work, &lwork,
+          iwork, &liwork, &info);
+      }
+
+      inline void heevd (
+        char const jobz, char const uplo, int const n,
+        traits::complex_f* a, int const lda,
+        float* w, traits::complex_f* work, int const lwork,
+        float* rwork, int const lrwork, int* iwork, int const liwork, int& info)
+      {
+        LAPACK_CHEEVD (
+          &jobz, &uplo, &n,
+          reinterpret_cast<fcomplex_t*>(a), &lda,
+          reinterpret_cast<fcomplex_t*>(w),
+          reinterpret_cast<fcomplex_t*>(work), &lwork,
+          rwork, &lrwork, iwork, &liwork, &info);
+      }
+
+      inline void heevd (
+        char const jobz, char const uplo, int const n,
+        traits::complex_d* a, int const lda,
+        double* w, traits::complex_d* work, int const lwork,
+        double* rwork, int const lrwork, int* iwork, int const liwork, int& info)
+      {
+        LAPACK_ZHEEVD (
+          &jobz, &uplo, &n,
+          reinterpret_cast<dcomplex_t*>(a), &lda,
+          reinterpret_cast<dcomplex_t*>(w),
+          reinterpret_cast<dcomplex_t*>(work), &lwork,
+          rwork, &lrwork, iwork, &liwork, &info);
+      }
+    } // namespace detail
+
+    namespace detail {
+
+      template <int N>
+      struct Heevd{};
+
+      /// Handling of workspace in the case of one workarray.
+      template <>
+      struct Heevd< 1 > {
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const uplo, int const n,
+          T* a, int const lda,
+          R* w, minimal_workspace, int& info) {
+
+          traits::detail::array<T> work( jobz=='N' ? 1+2*n : 1+6*n+2*n*n );
+          traits::detail::array<int> iwork( jobz=='N' ? 1 : 3+5*n );
+
+          heevd( jobz, uplo, n, a, lda, w,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (iwork), traits::vector_size (iwork),
+            info);
+        }
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const uplo, int const n,
+          T* a, int const lda,
+          R* w, optimal_workspace, int& info) {
+
+          traits::detail::array<int> iwork( jobz=='N' ? 1 : 3+5*n );
+
+          T workspace_query;
+          heevd( jobz, uplo, n, a, lda, w,
+            &workspace_query, -1,
+            traits::vector_storage (iwork), traits::vector_size (iwork),
+            info);
+
+          traits::detail::array<T> work( static_cast<int>( workspace_query ) );
+
+          heevd( jobz, uplo, n, a, lda, w,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (iwork), traits::vector_size (iwork),
+            info);
+        }
+        // Function that uses given workarrays
+        template <typename T, typename R, typename W, typename WI>
+        void operator() (
+          char const jobz, char const uplo, int const n,
+          T* a, int const lda,
+          R* w, std::pair<detail::workspace1<W>, detail::workspace1<WI> > work, int& info) {
+
+          assert (jobz=='N' ? 1+2*n : 1+6*n+2*n*n <= traits::vector_size (work.first.w_));
+          assert (jobz=='N' ? 1 : 3+5*n <= traits::vector_size (work.second.w_));
+
+          heevd( jobz, uplo, n, a, lda, w,
+            traits::vector_storage (work.first.w_), traits::vector_size (work.first.w_),
+            traits::vector_storage (work.second.w_), traits::vector_size (work.second.w_),
+            info);
+        }
+      };
+
+      /// Handling of workspace in the case of two workarrays.
+      template <>
+      struct Heevd< 2 > {
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const uplo, int const n,
+          T* a, int const lda,
+          R* w, minimal_workspace, int& info) {
+
+          traits::detail::array<T> work( jobz=='N' ? 1+n : 2*n+n*n );
+          traits::detail::array<R> rwork( jobz=='N' ? n : 1+5*n+2*n*n );
+          traits::detail::array<int> iwork( jobz=='N' ? 1 : 3+5*n );
+
+          heevd( jobz, uplo, n, a, lda, w,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (rwork), traits::vector_size (rwork),
+            traits::vector_storage (iwork), traits::vector_size (iwork),
+            info);
+        }
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const uplo, int const n,
+          T* a, int const lda,
+          R* w, optimal_workspace, int& info) {
+
+          traits::detail::array<R> rwork( jobz=='N' ? n : 1+5*n+2*n*n );
+          traits::detail::array<int> iwork( jobz=='N' ? 1 : 3+5*n );
+
+          T workspace_query;
+          heevd( jobz, uplo, n, a, lda, w,
+            &workspace_query, -1,
+            traits::vector_storage (rwork), traits::vector_size (rwork),
+            traits::vector_storage (iwork), traits::vector_size (iwork),
+            info);
+
+          traits::detail::array<T> work( static_cast<int>( traits::real( workspace_query ) ) );
+
+          heevd( jobz, uplo, n, a, lda, w,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (rwork), traits::vector_size (rwork),
+            traits::vector_storage (iwork), traits::vector_size (iwork),
+            info);
+        }
+        // Function that uses given workarrays
+        template <typename T, typename R, typename WC, typename WR, typename WI>
+        void operator() (
+          char const jobz, char const uplo, int const n,
+          T* a, int const lda,
+          R* w, std::pair<detail::workspace2<WC,WR>, detail::workspace1<WI> > work, int& info) {
+
+          assert (jobz=='N' ? 1+n : 2*n+n*n <= traits::vector_size (work.first.w_));
+          assert (jobz=='N' ? n : 1+5*n+2*n*n <= traits::vector_size (work.first.wr_));
+          assert (jobz=='N' ? 1 : 3+5*n <= traits::vector_size (work.second.w_));
+
+          heevd( jobz, uplo, n, a, lda, w,
+            traits::vector_storage (work.first.w_), traits::vector_size (work.first.w_),
+            traits::vector_storage (work.first.wr_), traits::vector_size (work.first.wr_),
+            traits::vector_storage (work.second.w_), traits::vector_size (work.second.w_),
+            info);
+        }
+      };
+    } // namespace detail
+
+    template <typename A, typename W, typename Work>
+    inline int heevd (
+      char jobz, char uplo, A& a,
+      W& w, Work work = optimal_workspace() ) {
+
+#ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
+      BOOST_STATIC_ASSERT((boost::is_same<
+        typename traits::matrix_traits<A>::matrix_structure, 
+        traits::general_t
+      >::value));
+#endif
+
+      int const n = traits::matrix_size1 (a);
+      assert (traits::matrix_size2 (a)==n);
+      assert (traits::vector_size (w)==n);
+      assert ( uplo=='U' || uplo=='L' );
+      assert ( jobz=='N' || jobz=='V' );
+
+      int info; 
+      detail::Heevd< n_workspace_args<typename A::value_type>::value >() (
+        jobz, uplo, n,
+        traits::matrix_storage (a),
+        traits::leading_dimension (a),
+        traits::vector_storage (w),
+        work,
+        info);
+      return info;
+    }
+  }
+
+}}}
+
+#endif
Added: sandbox/boost/numeric/bindings/lapack/heevx.hpp
==============================================================================
--- (empty file)
+++ sandbox/boost/numeric/bindings/lapack/heevx.hpp	2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
@@ -0,0 +1,321 @@
+/*
+ * 
+ * Copyright Toon Knapen, Karl Meerbergen & Kresimir Fresl 2003
+ * Copyright Thomas Klimpel 2008
+ *
+ * Distributed under the Boost Software License, Version 1.0.
+ * (See accompanying file LICENSE_1_0.txt or copy at
+ * http://www.boost.org/LICENSE_1_0.txt)
+ *
+ * KF acknowledges the support of the Faculty of Civil Engineering,
+ * University of Zagreb, Croatia.
+ *
+ */
+
+#ifndef BOOST_NUMERIC_BINDINGS_LAPACK_HEEVX_HPP
+#define BOOST_NUMERIC_BINDINGS_LAPACK_HEEVX_HPP
+
+#include <boost/numeric/bindings/traits/traits.hpp>
+#include <boost/numeric/bindings/traits/type.hpp>
+#include <boost/numeric/bindings/lapack/lapack.h>
+#include <boost/numeric/bindings/lapack/workspace.hpp>
+#include <boost/numeric/bindings/traits/detail/array.hpp>
+
+#ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
+#  include <boost/static_assert.hpp>
+#  include <boost/type_traits.hpp>
+#endif
+
+
+namespace boost { namespace numeric { namespace bindings {
+
+  namespace lapack {
+
+    ///////////////////////////////////////////////////////////////////
+    //
+    // Eigendecomposition of a complex Hermitian matrix A = Q * D * Q'
+    //
+    ///////////////////////////////////////////////////////////////////
+
+    /*
+     * heevx() computes selected eigenvalues and, optionally, eigenvectors
+     * of a complex Hermitian matrix A.  Eigenvalues and eigenvectors can
+     * be selected by specifying either a range of values or a range of
+     * indices for the desired eigenvalues.
+     *
+     * heevx() computes the eigendecomposition of a N x N matrix
+     * A = Q * D * Q',  where Q is a N x N unitary matrix and
+     * D is a diagonal matrix. The diagonal element D(i,i) is an
+     * eigenvalue of A and Q(:,i) is a corresponding eigenvector.
+     * The eigenvalues are stored in ascending order.
+     *
+     * On return of heevx, A is overwritten, z contains selected eigenvectors from Q
+     * and w contains selected eigenvalues from the main diagonal of D.
+     *
+     * int heevx (char jobz, char range, char uplo, A& a, T vl, T vu, int il, int iu, T abstol, int& m,
+     *            W& w, Z& z, IFail& ifail, Work work) {
+     *    jobz :  'V' : compute eigenvectors
+     *            'N' : do not compute eigenvectors
+     *    range : 'A': all eigenvalues will be found.
+     *            'V': all eigenvalues in the half-open interval (vl,vu] will be found.
+     *            'I': the il-th through iu-th eigenvalues will be found.
+     *    uplo :  'U' : only the upper triangular part of A is used on input.
+     *            'L' : only the lower triangular part of A is used on input.
+     */
+
+    namespace detail {
+
+      inline void heevx (
+        char const jobz, char const range, char const uplo, int const n,
+        float* a, int const lda,
+        float const vl, float const vu, int const il, int const iu,
+        float const abstol, int& m,
+        float* w, float* z, int const ldz, float* work, int const lwork,
+        int* iwork, int* ifail, int& info)
+      {
+        LAPACK_SSYEVX (
+          &jobz, &range, &uplo, &n,
+          a, &lda,
+          &vl, &vu, &il, &iu,
+          &abstol, &m,
+          w, z, &ldz, work, &lwork,
+          iwork, ifail, &info);
+      }
+
+      inline void heevx (
+        char const jobz, char const range, char const uplo, int const n,
+        double* a, int const lda,
+        double const vl, double const vu, int const il, int const iu,
+        double const abstol, int& m,
+        double* w, double* z, int const ldz, double* work, int const lwork,
+        int* iwork, int* ifail, int& info)
+      {
+        LAPACK_DSYEVX (
+          &jobz, &range, &uplo, &n,
+          a, &lda,
+          &vl, &vu, &il, &iu,
+          &abstol, &m,
+          w, z, &ldz, work, &lwork,
+          iwork, ifail, &info);
+      }
+
+      inline void heevx (
+        char const jobz, char const range, char const uplo, int const n,
+        traits::complex_f* a, int const lda,
+        float const vl, float const vu, int const il, int const iu,
+        float const abstol, int& m,
+        float* w, traits::complex_f* z, int const ldz, traits::complex_f* work, int const lwork,
+        float* rwork, int* iwork, int* ifail, int& info)
+      {
+        LAPACK_CHEEVX (
+          &jobz, &range, &uplo, &n,
+          reinterpret_cast<fcomplex_t*>(a), &lda,
+          &vl, &vu, &il, &iu, &abstol, &m, w,
+          reinterpret_cast<fcomplex_t*>(z), &ldz,
+          reinterpret_cast<fcomplex_t*>(work), &lwork,
+          rwork, iwork, ifail, &info);
+      }
+
+      inline void heevx (
+        char const jobz, char const range, char const uplo, int const n,
+        traits::complex_d* a, int const lda,
+        double const vl, double const vu, int const il, int const iu,
+        double const abstol, int& m,
+        double* w, traits::complex_d* z, int const ldz, traits::complex_d* work, int const lwork,
+        double* rwork, int* iwork, int* ifail, int& info)
+      {
+        LAPACK_ZHEEVX (
+          &jobz, &range, &uplo, &n,
+          reinterpret_cast<dcomplex_t*>(a), &lda,
+          &vl, &vu, &il, &iu, &abstol, &m, w,
+          reinterpret_cast<dcomplex_t*>(z), &ldz,
+          reinterpret_cast<dcomplex_t*>(work), &lwork,
+          rwork, iwork, ifail, &info);
+      }
+    } // namespace detail
+
+    namespace detail {
+
+      template <int N>
+      struct Heevx{};
+
+      /// Handling of workspace in the case of one workarray.
+      template <>
+      struct Heevx< 1 > {
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const range, char const uplo, int const n,
+          T* a, int const lda,
+          R vl, R vu, int const il, int const iu,
+          R abstol, int& m,
+          R* w, T* z, int const ldz, minimal_workspace, int* ifail, int& info) {
+
+          traits::detail::array<T> work( 8*n );
+          traits::detail::array<int> iwork( 5*n );
+
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (iwork),
+            ifail, info);
+        }
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const range, char const uplo, int const n,
+          T* a, int const lda,
+          R vl, R vu, int const il, int const iu,
+          R abstol, int& m,
+          R* w, T* z, int const ldz, optimal_workspace, int* ifail, int& info) {
+
+          traits::detail::array<int> iwork( 5*n );
+
+          T workspace_query;
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            &workspace_query, -1,
+            traits::vector_storage (iwork),
+            ifail, info);
+
+          traits::detail::array<T> work( static_cast<int>( workspace_query ) );
+
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (iwork),
+            ifail, info);
+        }
+        // Function that uses given workarrays
+        template <typename T, typename R, typename W, typename WI>
+        void operator() (
+          char const jobz, char const range, char const uplo, int const n,
+          T* a, int const lda,
+          R vl, R vu, int const il, int const iu,
+          R abstol, int& m,
+          R* w, T* z, int const ldz, std::pair<detail::workspace1<W>, detail::workspace1<WI> > work, int* ifail, int& info) {
+
+          assert (8*n <= traits::vector_size (work.first.w_));
+          assert (5*n <= traits::vector_size (work.second.w_));
+
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            traits::vector_storage (work.first.w_), traits::vector_size (work.first.w_),
+            traits::vector_storage (work.second.w_),
+            ifail, info);
+        }
+      };
+
+      /// Handling of workspace in the case of two workarrays.
+      template <>
+      struct Heevx< 2 > {
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const range, char const uplo, int const n,
+          T* a, int const lda,
+          R vl, R vu, int const il, int const iu,
+          R abstol, int& m,
+          R* w, T* z, int const ldz, minimal_workspace, int* ifail, int& info) {
+
+          traits::detail::array<T> work( 2*n );
+          traits::detail::array<R> rwork( 7*n );
+          traits::detail::array<int> iwork( 5*n );
+
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (rwork),
+            traits::vector_storage (iwork),
+            ifail, info);
+        }
+        // Function that allocates temporary arrays
+        template <typename T, typename R>
+        void operator() (
+          char const jobz, char const range, char const uplo, int const n,
+          T* a, int const lda,
+          R vl, R vu, int const il, int const iu,
+          R abstol, int& m,
+          R* w, T* z, int const ldz, optimal_workspace, int* ifail, int& info) {
+
+          traits::detail::array<R> rwork( 7*n );
+          traits::detail::array<int> iwork( 5*n );
+
+          T workspace_query;
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            &workspace_query, -1,
+            traits::vector_storage (rwork),
+            traits::vector_storage (iwork),
+            ifail, info);
+
+          traits::detail::array<T> work( static_cast<int>( traits::real( workspace_query ) ) );
+
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            traits::vector_storage (work), traits::vector_size (work),
+            traits::vector_storage (rwork),
+            traits::vector_storage (iwork),
+            ifail, info);
+        }
+        // Function that uses given workarrays
+        template <typename T, typename R, typename WC, typename WR, typename WI>
+        void operator() (
+          char const jobz, char const range, char const uplo, int const n,
+          T* a, int const lda,
+          R vl, R vu, int const il, int const iu,
+          R abstol, int& m,
+          R* w, T* z, int const ldz, std::pair<detail::workspace2<WC,WR>, detail::workspace1<WI> > work, int* ifail, int& info) {
+
+          assert (2*n <= traits::vector_size (work.first.w_));
+          assert (7*n <= traits::vector_size (work.first.wr_));
+          assert (5*n <= traits::vector_size (work.second.w_));
+
+          heevx( jobz, range, uplo, n, a, lda, vl, vu, il, iu, abstol, m, w, z, ldz,
+            traits::vector_storage (work.first.w_), traits::vector_size (work.first.w_),
+            traits::vector_storage (work.first.wr_),
+            traits::vector_storage (work.second.w_),
+            ifail, info);
+        }
+      };
+    } // namespace detail
+
+    template <typename A, typename T, typename W, typename Z, typename IFail, typename Work>
+    inline int heevx (
+      char jobz, char range, A& a, T vl, T vu, int il, int iu, T abstol, int& m,
+      //char jobz, char range, char uplo, A& a, T vl, T vu, int il, int iu, T abstol, int& m,
+      W& w, Z& z, IFail& ifail, Work work = optimal_workspace() ) {
+
+#ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
+      typedef typename A::value_type                               value_type ;
+      typedef typename traits::type_traits< value_type >::real_type real_type ;
+      BOOST_STATIC_ASSERT((boost::is_same<
+        typename traits::matrix_traits<A>::matrix_structure, 
+        traits::hermitian_t
+      >::value || (boost::is_same<
+        typename traits::matrix_traits<A>::matrix_structure,
+        traits::symmetric_t
+      >::value && boost::is_same<value_type, real_type>::value)));
+#endif
+
+      int const n = traits::matrix_size1 (a);
+      assert (traits::matrix_size2 (a)==n);
+      assert (traits::vector_size (w)==n);
+      assert (n == traits::vector_size (ifail));
+      assert ( range=='A' || range=='V' || range=='I' );
+      char uplo = traits::matrix_uplo_tag (a);
+      assert ( uplo=='U' || uplo=='L' );
+      assert ( jobz=='N' || jobz=='V' );
+
+      int info; 
+      detail::Heevx< n_workspace_args<typename A::value_type>::value >() (
+        jobz, range, uplo, n,
+        traits::matrix_storage (a),
+        traits::leading_dimension (a),
+        vl, vu, il, iu, abstol, m,
+        traits::vector_storage (w),
+        traits::matrix_storage (z),
+        traits::leading_dimension (z),
+        work,
+        traits::vector_storage (ifail),
+        info);
+      return info;
+    }
+  }
+
+}}}
+
+#endif
Modified: sandbox/boost/numeric/bindings/lapack/lapack.h
==============================================================================
--- sandbox/boost/numeric/bindings/lapack/lapack.h	(original)
+++ sandbox/boost/numeric/bindings/lapack/lapack.h	2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
@@ -323,6 +323,52 @@
                      dcomplex_t* work, const int * lwork, double* rwork,
                      int* info );
 
+  
+  void LAPACK_SSYEVD( const char* jobz, const char* uplo, const int* n,
+                      float* a, const int* lda, float* w,
+                      float* work, const int* lwork,
+                      int* iwork, const int* liwork, int* info);
+
+  void LAPACK_DSYEVD( const char* jobz, const char* uplo, const int* n,
+                      double* a, const int* lda, double* w,
+                      double* work, const int* lwork,
+                      int* iwork, const int* liwork, int* info);
+
+  void LAPACK_CHEEVD( const char* jobz, const char* uplo, const int* n,
+                      fcomplex_t* a, const int* lda, float* w,
+                      fcomplex_t* work, const int* lwork, float* rwork, const int* lrwork,
+                      int* iwork, const int* liwork, int* info);
+
+  void LAPACK_ZHEEVD( const char* jobz, const char* uplo, const int* n,
+                      dcomplex_t* a, const int* lda, double* w,
+                      dcomplex_t* work, const int* lwork, double* rwork, const int* lrwork,
+                      int* iwork, const int* liwork, int* info);
+
+
+  void LAPACK_SSYEVX( const char* jobz, const char* range, const char* uplo, const int* n,
+                      float* a, const int* lda, const float* vl, const float* vu, const int* il, const int* iu,
+                      const float* abstol, int* m, float* w, float* z, const int* ldz,
+                      float* work, const int* lwork,
+                      int* iwork, int* ifail, int* info);
+
+  void LAPACK_DSYEVX( const char* jobz, const char* range, const char* uplo, const int* n,
+                      double* a, const int* lda, const double* vl, const double* vu, const int* il, const int* iu,
+                      const double* abstol, int* m, double* w, double* z, const int* ldz,
+                      double* work, const int* lwork,
+                      int* iwork, int* ifail, int* info);
+
+  void LAPACK_CHEEVX( const char* jobz, const char* range, const char* uplo, const int* n,
+                      fcomplex_t* a, const int* lda, const float* vl, const float* vu, const int* il, const int* iu,
+                      const float* abstol, int* m, float* w, fcomplex_t* z, const int* ldz,
+                      fcomplex_t* work, const int* lwork, float* rwork,
+                      int* iwork, int* ifail, int* info);
+
+  void LAPACK_ZHEEVX( const char* jobz, const char* range, const char* uplo, const int* n,
+                      dcomplex_t* a, const int* lda, const double* vl, const double* vu, const int* il, const int* iu,
+                      const double* abstol, int* m, double* w, dcomplex_t* z, const int* ldz,
+                      dcomplex_t* work, const int* lwork, double* rwork,
+                      int* iwork, int* ifail, int* info);
+
 
   void LAPACK_CTREVC( const char* side, const char* howmny, const logical_t* select, const int *n,
                      fcomplex_t* t, const int * ldt, fcomplex_t* vl, const int* ldvl,
Modified: sandbox/boost/numeric/bindings/lapack/lapack_names.h
==============================================================================
--- sandbox/boost/numeric/bindings/lapack/lapack_names.h	(original)
+++ sandbox/boost/numeric/bindings/lapack/lapack_names.h	2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
@@ -157,6 +157,16 @@
 #define LAPACK_CHEEV FORTRAN_ID( cheev )
 #define LAPACK_ZHEEV FORTRAN_ID( zheev )
 
+#define LAPACK_SSYEVD FORTRAN_ID( ssyevd )
+#define LAPACK_DSYEVD FORTRAN_ID( dsyevd )
+#define LAPACK_CHEEVD FORTRAN_ID( cheevd )
+#define LAPACK_ZHEEVD FORTRAN_ID( zheevd )
+
+#define LAPACK_SSYEVX FORTRAN_ID( ssyevx )
+#define LAPACK_DSYEVX FORTRAN_ID( dsyevx )
+#define LAPACK_CHEEVX FORTRAN_ID( cheevx )
+#define LAPACK_ZHEEVX FORTRAN_ID( zheevx )
+
 
 #define LAPACK_STREVC FORTRAN_ID( strevc )
 #define LAPACK_DTREVC FORTRAN_ID( dtrevc )
Added: sandbox/boost/numeric/bindings/lapack/syevd.hpp
==============================================================================
--- (empty file)
+++ sandbox/boost/numeric/bindings/lapack/syevd.hpp	2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
@@ -0,0 +1,33 @@
+// 
+// Copyright (c) Thomas Klimpel 2008
+//
+// Distributed under the Boost Software License, Version 1.0.
+// (See accompanying file LICENSE_1_0.txt or copy at
+// http://www.boost.org/LICENSE_1_0.txt)
+//
+
+#ifndef BOOST_NUMERIC_BINDINGS_LAPACK_SYEVD_HPP
+#define BOOST_NUMERIC_BINDINGS_LAPACK_SYEVD_HPP
+
+#include <boost/numeric/bindings/lapack/heevd.hpp>
+
+namespace boost { namespace numeric { namespace bindings { 
+  namespace lapack {
+
+    template <typename A, typename W, typename Work>
+    inline int syevd (
+      char jobz, char uplo, A& a,
+      W& w, Work work = optimal_workspace() ) {
+
+#ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
+      typedef typename A::value_type                               value_type ;
+      typedef typename traits::type_traits< value_type >::real_type real_type ;
+      BOOST_STATIC_ASSERT((boost::is_same<value_type, real_type>::value));
+#endif
+
+      return heevd (jobz, uplo, a, w, work);
+    }
+  }
+}}}
+
+#endif
Added: sandbox/boost/numeric/bindings/lapack/syevx.hpp
==============================================================================
--- (empty file)
+++ sandbox/boost/numeric/bindings/lapack/syevx.hpp	2008-05-15 15:00:57 EDT (Thu, 15 May 2008)
@@ -0,0 +1,32 @@
+// 
+// Copyright (c) Thomas Klimpel 2008
+//
+// Distributed under the Boost Software License, Version 1.0.
+// (See accompanying file LICENSE_1_0.txt or copy at
+// http://www.boost.org/LICENSE_1_0.txt)
+//
+
+#ifndef BOOST_NUMERIC_BINDINGS_LAPACK_SYEVX_HPP
+#define BOOST_NUMERIC_BINDINGS_LAPACK_SYEVX_HPP
+
+#include <boost/numeric/bindings/lapack/heevx.hpp>
+
+namespace boost { namespace numeric { namespace bindings { 
+  namespace lapack {
+    template <typename A, typename T, typename W, typename Z, typename IFail, typename Work>
+    int syevx (
+      char jobz, char range, A& a, T vl, T vu, int il, int iu, T abstol, int& m,
+      W& w, Z& z, IFail& ifail, Work work = optimal_workspace() ) {
+
+#ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
+      typedef typename A::value_type                               value_type ;
+      typedef typename traits::type_traits< value_type >::real_type real_type ;
+      BOOST_STATIC_ASSERT((boost::is_same<value_type, real_type>::value));
+#endif
+
+      return heevx (jobz, range, a, vl, vu, il, iu, abstol, m, w, z, ifail, work);
+    }
+  }
+}}}
+
+#endif