clapack.hpp

Go to the documentation of this file.
00001 /*
00002  * 
00003  * Copyright (c) Kresimir Fresl 2002 
00004  *
00005  * Distributed under the Boost Software License, Version 1.0.
00006  * (See accompanying file LICENSE_1_0.txt or copy at
00007  * http://www.boost.org/LICENSE_1_0.txt)
00008  *
00009  * Author acknowledges the support of the Faculty of Civil Engineering, 
00010  * University of Zagreb, Croatia.
00011  *
00012  */
00013 
00014 #ifndef BOOST_NUMERIC_BINDINGS_CLAPACK_HPP
00015 #define BOOST_NUMERIC_BINDINGS_CLAPACK_HPP
00016 
00017 #include <cassert>
00018 #include <new>
00019 
00020 #include <boost/numeric/bindings/traits/traits.hpp>
00021 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00022 #  include <boost/numeric/bindings/traits/detail/symm_herm_traits.hpp>
00023 #endif 
00024 #include <boost/numeric/bindings/atlas/cblas_enum.hpp>
00025 
00026 // see libs/numeric/bindings/atlas/doc/index.html, section 2.5.2
00027 //#define BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00028 
00029 #include <boost/numeric/bindings/atlas/clapack_overloads.hpp>
00030 
00031 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00032 #  include <boost/static_assert.hpp>
00033 #  include <boost/type_traits/same_traits.hpp>
00034 #endif
00035 
00036 
00037 namespace boost { namespace numeric { namespace bindings { 
00038 
00039   namespace atlas {
00040 
00042     //
00043     // general system of linear equations A * X = B
00044     // 
00046 
00047     // gesv(): 'driver' function 
00048     //
00049     // [comments from 'clapack_dgesv.c':]
00050     /* clapack_xgesv computes the solution to a system of linear equations
00051      *   A * X = B,
00052      * where A is an N-by-N matrix and X and B are N-by-NRHS matrices.
00053      */
00054     // [but ATLAS FAQ says:]
00055     /* What's the deal with the RHS in the row-major factorization/solves?
00056      * Most users are confused by the row major factorization and related 
00057      * solves. The right-hand side vectors are probably the biggest source 
00058      * of confusion. The RHS array does not represent a matrix in the 
00059      * mathematical sense, it is instead a pasting together of the various 
00060      * RHS into one array for calling convenience. As such, RHS vectors are 
00061      * always stored contiguously, regardless of the row/col major that is 
00062      * chosen. This means that ldb/ldx is always independent of NRHS, and 
00063      * dependant on N, regardless of the row/col major setting. 
00064      */ 
00065     // That is, it seems that, if B is row-major, it should be NRHS-by-N, 
00066     // and RHS vectors should be its rows, not columns. 
00067     //
00068     // [comments from 'clapack_dgesv.c':]
00069     /* The LU factorization used to factor A is dependent on the Order 
00070      * parameter, as detailed in the leading comments of clapack_dgetrf.
00071      * The factored form of A is then used to solve the system of equations 
00072      *   A * X = B.
00073      * A is overwritten with the appropriate LU factorization, and B [...]
00074      * is overwritten with the solution X on output.
00075      */
00076     // If B is row-major, solution vectors are its rows. 
00077     template <typename MatrA, typename MatrB, typename IVec>
00078     inline
00079     int gesv (MatrA& a, IVec& ipiv, MatrB& b) {
00080 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00081       BOOST_STATIC_ASSERT((boost::is_same<
00082         typename traits::matrix_traits<MatrA>::matrix_structure, 
00083         traits::general_t
00084       >::value)); 
00085       BOOST_STATIC_ASSERT((boost::is_same<
00086         typename traits::matrix_traits<MatrB>::matrix_structure, 
00087         traits::general_t
00088       >::value)); 
00089 
00090       BOOST_STATIC_ASSERT((boost::is_same<
00091         typename traits::matrix_traits<MatrA>::ordering_type,
00092         typename traits::matrix_traits<MatrB>::ordering_type
00093       >::value)); 
00094 #endif 
00095 
00096       CBLAS_ORDER const stor_ord
00097         = enum_cast<CBLAS_ORDER const>
00098         (storage_order<
00099 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00100            typename traits::matrix_traits<MatrA>::ordering_type
00101 #else
00102            typename MatrA::orientation_category 
00103 #endif 
00104          >::value); 
00105 
00106       int const n = traits::matrix_size1 (a);
00107       int const nrhs = stor_ord == CblasColMajor
00108         ? traits::matrix_size2 (b)
00109         : traits::matrix_size1 (b); 
00110       assert (n == traits::matrix_size2 (a)); 
00111       assert (n == (stor_ord == CblasColMajor
00112                     ? traits::matrix_size1 (b)
00113                     : traits::matrix_size2 (b))); 
00114       assert (n == traits::vector_size (ipiv)); 
00115 
00116       return detail::gesv (stor_ord, n, nrhs, 
00117                            traits::matrix_storage (a), 
00118                            traits::leading_dimension (a),
00119                            traits::vector_storage (ipiv),  
00120                            traits::matrix_storage (b),
00121                            traits::leading_dimension (b));
00122     }
00123 
00124     template <typename MatrA, typename MatrB>
00125     inline
00126     int gesv (MatrA& a, MatrB& b) {
00127       // with 'internal' pivot vector
00128 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00129       BOOST_STATIC_ASSERT((boost::is_same<
00130         typename traits::matrix_traits<MatrA>::matrix_structure, 
00131         traits::general_t
00132       >::value)); 
00133       BOOST_STATIC_ASSERT((boost::is_same<
00134         typename traits::matrix_traits<MatrB>::matrix_structure, 
00135         traits::general_t
00136       >::value)); 
00137 
00138       BOOST_STATIC_ASSERT((boost::is_same<
00139         typename traits::matrix_traits<MatrA>::ordering_type,
00140         typename traits::matrix_traits<MatrB>::ordering_type
00141       >::value)); 
00142 #endif 
00143 
00144       CBLAS_ORDER const stor_ord
00145         = enum_cast<CBLAS_ORDER const>
00146         (storage_order<
00147 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00148            typename traits::matrix_traits<MatrA>::ordering_type
00149 #else
00150            typename MatrA::orientation_category 
00151 #endif 
00152          >::value); 
00153 
00154       int const n = traits::matrix_size1 (a);
00155       int const nrhs = stor_ord == CblasColMajor
00156         ? traits::matrix_size2 (b)
00157         : traits::matrix_size1 (b); 
00158       assert (n == traits::matrix_size2 (a)); 
00159       assert (n == (stor_ord == CblasColMajor
00160                     ? traits::matrix_size1 (b)
00161                     : traits::matrix_size2 (b))); 
00162 
00163       int *ipiv = new (std::nothrow) int[n]; 
00164       int ierr = -101;  
00165       // clapack_dgesv() errors: 
00166       //   if (ierr == 0), successful
00167       //   if (ierr < 0), the -ierr argument had an illegal value
00168       //   -- we will use -101 if allocation fails
00169       //   if (ierr > 0), U(i-1,i-1) (or L(i-1,i-1)) is exactly zero 
00170  
00171       if (ipiv) {
00172         ierr = detail::gesv (stor_ord, n, nrhs, 
00173                              traits::matrix_storage (a), 
00174                              traits::leading_dimension (a),
00175                              ipiv,  
00176                              traits::matrix_storage (b),
00177                              traits::leading_dimension (b));
00178         delete[] ipiv; 
00179       }
00180       return ierr; 
00181     }
00182 
00183     template <typename MatrA, typename MatrB>
00184     inline
00185     int lu_solve (MatrA& a, MatrB& b) {
00186       return gesv (a, b); 
00187     }
00188 
00189 
00190     // getrf(): LU factorization of A
00191     // [comments from 'clapack_dgetrf.c':]
00192     /* Computes one of two LU factorizations based on the setting of 
00193      * the Order parameter, as follows:
00194      * ---------------------------------------------------------------
00195      *                     Order == CblasColMajor
00196      * Column-major factorization of form
00197      *   A = P * L * U
00198      * where P is a row-permutation matrix, L is lower triangular with
00199      * unit diagonal elements (lower trapezoidal if M > N), and U is 
00200      * upper triangular (upper trapezoidal if M < N).
00201      * ---------------------------------------------------------------
00202      *                     Order == CblasRowMajor
00203      * Row-major factorization of form
00204      *   A = P * L * U
00205      * where P is a column-permutation matrix, L is lower triangular 
00206      * (lower trapezoidal if M > N), and U is upper triangular with 
00207      * unit diagonals (upper trapezoidal if M < N).
00208      */
00209     template <typename MatrA, typename IVec> 
00210     inline
00211     int getrf (MatrA& a, IVec& ipiv) {
00212 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00213       BOOST_STATIC_ASSERT((boost::is_same<
00214         typename traits::matrix_traits<MatrA>::matrix_structure, 
00215         traits::general_t
00216       >::value)); 
00217 #endif 
00218 
00219       CBLAS_ORDER const stor_ord
00220         = enum_cast<CBLAS_ORDER const>
00221         (storage_order<
00222 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00223            typename traits::matrix_traits<MatrA>::ordering_type
00224 #else
00225            typename MatrA::orientation_category 
00226 #endif 
00227          >::value); 
00228 
00229       int const m = traits::matrix_size1 (a);
00230       int const n = traits::matrix_size2 (a); 
00231       assert (traits::vector_size (ipiv) == (m < n ? m : n)); 
00232 
00233       return detail::getrf (stor_ord, m, n, 
00234                             traits::matrix_storage (a), 
00235                             traits::leading_dimension (a),
00236                             traits::vector_storage (ipiv)); 
00237     }
00238 
00239     template <typename MatrA, typename IVec> 
00240     inline
00241     int lu_factor (MatrA& a, IVec& ipiv) {
00242       return getrf (a, ipiv); 
00243     }
00244 
00245 
00246     // getrs(): solves a system of linear equations
00247     //          A * X = B  or  A' * X = B
00248     //          using the LU factorization previously computed by getrf()
00249     template <typename MatrA, typename MatrB, typename IVec>
00250     inline
00251     int getrs (CBLAS_TRANSPOSE const Trans, 
00252                MatrA const& a, IVec const& ipiv, MatrB& b) 
00253     {
00254 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00255       BOOST_STATIC_ASSERT((boost::is_same<
00256         typename traits::matrix_traits<MatrA>::matrix_structure, 
00257         traits::general_t
00258       >::value)); 
00259       BOOST_STATIC_ASSERT((boost::is_same<
00260         typename traits::matrix_traits<MatrB>::matrix_structure, 
00261         traits::general_t
00262       >::value)); 
00263 
00264       BOOST_STATIC_ASSERT((boost::is_same<
00265         typename traits::matrix_traits<MatrA>::ordering_type,
00266         typename traits::matrix_traits<MatrB>::ordering_type
00267       >::value)); 
00268 #endif 
00269 
00270       assert (Trans == CblasNoTrans 
00271               || Trans == CblasTrans 
00272               || Trans == CblasConjTrans); 
00273 
00274       CBLAS_ORDER const stor_ord
00275         = enum_cast<CBLAS_ORDER const>
00276         (storage_order<
00277 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00278            typename traits::matrix_traits<MatrA>::ordering_type
00279 #else
00280            typename MatrA::orientation_category 
00281 #endif 
00282          >::value); 
00283 
00284       int const n = traits::matrix_size1 (a);
00285       int const nrhs = stor_ord == CblasColMajor
00286         ? traits::matrix_size2 (b)
00287         : traits::matrix_size1 (b); 
00288       assert (n == traits::matrix_size2 (a)); 
00289       assert (n == (stor_ord == CblasColMajor
00290                     ? traits::matrix_size1 (b)
00291                     : traits::matrix_size2 (b))); 
00292       assert (n == traits::vector_size (ipiv)); 
00293       
00294       return detail::getrs (stor_ord, Trans, n, nrhs, 
00295 #ifndef BOOST_NO_FUNCTION_TEMPLATE_ORDERING
00296                             traits::matrix_storage (a), 
00297 #else
00298                             traits::matrix_storage_const (a), 
00299 #endif 
00300                             traits::leading_dimension (a),
00301 #ifndef BOOST_NO_FUNCTION_TEMPLATE_ORDERING
00302                             traits::vector_storage (ipiv),  
00303 #else
00304                             traits::vector_storage_const (ipiv),  
00305 #endif
00306                             traits::matrix_storage (b),
00307                             traits::leading_dimension (b)); 
00308     }
00309 
00310     // getrs(): solves A * X = B (after getrf())
00311     template <typename MatrA, typename MatrB, typename IVec>
00312     inline
00313     int getrs (MatrA const& a, IVec const& ipiv, MatrB& b) {
00314       return getrs (CblasNoTrans, a, ipiv, b); 
00315     }
00316 
00317     template <typename MatrA, typename MatrB, typename IVec>
00318     inline
00319     int lu_substitute (MatrA const& a, IVec const& ipiv, MatrB& b) {
00320       return getrs (CblasNoTrans, a, ipiv, b); 
00321     }
00322 
00323 
00324     // getri(): computes the inverse of a matrix A 
00325     //          using the LU factorization previously computed by getrf() 
00326     template <typename MatrA, typename IVec> 
00327     inline
00328     int getri (MatrA& a, IVec const& ipiv) {
00329 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00330       BOOST_STATIC_ASSERT((boost::is_same<
00331         typename traits::matrix_traits<MatrA>::matrix_structure, 
00332         traits::general_t
00333       >::value)); 
00334 #endif 
00335 
00336       CBLAS_ORDER const stor_ord
00337         = enum_cast<CBLAS_ORDER const>
00338         (storage_order<
00339 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00340            typename traits::matrix_traits<MatrA>::ordering_type
00341 #else
00342            typename MatrA::orientation_category 
00343 #endif 
00344          >::value); 
00345 
00346       int const n = traits::matrix_size1 (a);
00347       assert (n == traits::matrix_size2 (a)); 
00348       assert (traits::vector_size (ipiv) == n); 
00349 
00350       return detail::getri (stor_ord, n, 
00351                             traits::matrix_storage (a), 
00352                             traits::leading_dimension (a),
00353 #ifndef BOOST_NO_FUNCTION_TEMPLATE_ORDERING
00354                             traits::vector_storage (ipiv)
00355 #else
00356                             traits::vector_storage_const (ipiv)
00357 #endif
00358                             ); 
00359     }
00360 
00361     template <typename MatrA, typename IVec> 
00362     inline
00363     int lu_invert (MatrA& a, IVec& ipiv) {
00364       return getri (a, ipiv); 
00365     }
00366 
00367 
00368 
00370     //
00371     // system of linear equations A * X = B
00372     // with A symmetric or Hermitian positive definite matrix
00373     //
00375 
00376 #ifndef BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00377     // posv(): 'driver' function 
00378     //
00379     // [from 'dposv.f' (slightly edited):]
00380     /* XPOSV computes the solution to a system of linear equations
00381      *    A * X = B,
00382      * where A is an N-by-N symmetric/Hermitian positive definite matrix 
00383      * and X and B are N-by-NRHS matrices. [See also comments of gesv().]
00384      *
00385      * A -- On entry, the symmetric/Hermitian matrix A.  
00386      * If UPLO = 'U', the leading N-by-N upper triangular part of A 
00387      * contains the upper triangular part of the matrix A, and the 
00388      * strictly lower triangular part of A is not referenced.  
00389      * If UPLO = 'L', the leading N-by-N lower triangular part of A 
00390      * contains the lower triangular part of the matrix A, and the 
00391      * strictly upper triangular part of A is not referenced.
00392      *
00393      * On exit, if INFO = 0, the factor U or L from the Cholesky
00394      * factorization A = U**T*U or A = L*L**T
00395      * [or A = U**H*U or A = L*L**H]. 
00396      *
00397      * B -- On entry, the right hand side matrix B.
00398      * On exit, if INFO = 0, the solution matrix X.
00399      */
00400     namespace detail {
00401 
00402       template <typename SymmA, typename MatrB>
00403       inline
00404       int posv (CBLAS_UPLO const uplo, SymmA& a, MatrB& b) {
00405 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00406         BOOST_STATIC_ASSERT((boost::is_same<
00407           typename traits::matrix_traits<SymmA>::ordering_type,
00408           typename traits::matrix_traits<MatrB>::ordering_type
00409         >::value)); 
00410 #endif 
00411 
00412         CBLAS_ORDER const stor_ord
00413           = enum_cast<CBLAS_ORDER const>
00414           (storage_order<
00415 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00416             typename traits::matrix_traits<SymmA>::ordering_type
00417 #else
00418             typename SymmA::orientation_category 
00419 #endif 
00420            >::value); 
00421 
00422         int const n = traits::matrix_size1 (a);
00423         int const nrhs = stor_ord == CblasColMajor
00424           ? traits::matrix_size2 (b)
00425           : traits::matrix_size1 (b); 
00426         assert (n == traits::matrix_size2 (a)); 
00427         assert (n == (stor_ord == CblasColMajor
00428                       ? traits::matrix_size1 (b)
00429                       : traits::matrix_size2 (b))); 
00430 
00431         return posv (stor_ord, uplo, n, nrhs, 
00432                      traits::matrix_storage (a), 
00433                      traits::leading_dimension (a),
00434                      traits::matrix_storage (b),
00435                      traits::leading_dimension (b));
00436       }
00437 
00438     } // detail 
00439 
00440     template <typename SymmA, typename MatrB>
00441     inline
00442     int posv (CBLAS_UPLO const uplo, SymmA& a, MatrB& b) {
00443 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
00444       BOOST_STATIC_ASSERT((boost::is_same<
00445         typename traits::matrix_traits<SymmA>::matrix_structure, 
00446         traits::general_t
00447       >::value));
00448       BOOST_STATIC_ASSERT((boost::is_same<
00449         typename traits::matrix_traits<MatrB>::matrix_structure, 
00450         traits::general_t
00451       >::value));
00452 #endif
00453       assert (uplo == CblasUpper || uplo == CblasLower); 
00454       return detail::posv (uplo, a, b); 
00455     }
00456 
00457     template <typename SymmA, typename MatrB>
00458     inline
00459     int posv (SymmA& a, MatrB& b) {
00460 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00461       typedef typename traits::matrix_traits<SymmA>::value_type val_t;
00462       BOOST_STATIC_ASSERT((boost::is_same<
00463         typename traits::matrix_traits<SymmA>::matrix_structure, 
00464         typename traits::detail::symm_herm_t<val_t>::type
00465       >::value)); 
00466       BOOST_STATIC_ASSERT((boost::is_same<
00467         typename traits::matrix_traits<MatrB>::matrix_structure, 
00468         traits::general_t
00469       >::value)); 
00470 #endif 
00471 
00472       CBLAS_UPLO const uplo
00473         = enum_cast<CBLAS_UPLO const>
00474         (uplo_triang<
00475 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00476           typename traits::matrix_traits<SymmA>::uplo_type
00477 #else
00478           typename SymmA::packed_category 
00479 #endif 
00480          >::value); 
00481       
00482       return detail::posv (uplo, a, b); 
00483     }
00484 
00485     template <typename SymmA, typename MatrB>
00486     inline
00487     int cholesky_solve (SymmA& a, MatrB& b) { return posv (a, b); }
00488 #endif // BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00489 
00490 
00491     // potrf(): Cholesky factorization of A 
00492     namespace detail {
00493 
00494       template <typename SymmA>
00495       inline
00496       int potrf (CBLAS_UPLO const uplo, SymmA& a) {
00497 
00498         CBLAS_ORDER const stor_ord
00499           = enum_cast<CBLAS_ORDER const>
00500           (storage_order<
00501 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00502             typename traits::matrix_traits<SymmA>::ordering_type
00503 #else
00504             typename SymmA::orientation_category 
00505 #endif 
00506            >::value); 
00507 
00508         int const n = traits::matrix_size1 (a);
00509         assert (n == traits::matrix_size2 (a)); 
00510 
00511         return potrf (stor_ord, uplo, n, 
00512                       traits::matrix_storage (a), 
00513                       traits::leading_dimension (a));
00514       }
00515 
00516     } // detail 
00517 
00518     template <typename SymmA>
00519     inline
00520     int potrf (CBLAS_UPLO const uplo, SymmA& a) {
00521 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
00522       BOOST_STATIC_ASSERT((boost::is_same<
00523         typename traits::matrix_traits<SymmA>::matrix_structure, 
00524         traits::general_t
00525       >::value));
00526 #endif
00527       assert (uplo == CblasUpper || uplo == CblasLower); 
00528       return detail::potrf (uplo, a); 
00529     } 
00530 
00531     template <typename SymmA>
00532     inline
00533     int potrf (SymmA& a) {
00534 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00535       typedef typename traits::matrix_traits<SymmA>::value_type val_t;
00536       BOOST_STATIC_ASSERT((boost::is_same<
00537         typename traits::matrix_traits<SymmA>::matrix_structure, 
00538         typename traits::detail::symm_herm_t<val_t>::type
00539       >::value)); 
00540 #endif 
00541 
00542       CBLAS_UPLO const uplo
00543         = enum_cast<CBLAS_UPLO const>
00544         (uplo_triang<
00545 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00546           typename traits::matrix_traits<SymmA>::uplo_type
00547 #else
00548           typename SymmA::packed_category 
00549 #endif 
00550          >::value); 
00551       
00552       return detail::potrf (uplo, a); 
00553     }
00554 
00555     template <typename SymmA>
00556     inline
00557     int cholesky_factor (SymmA& a) { return potrf (a); }
00558 
00559 
00560     // potrs(): solves a system of linear equations A * X = B
00561     //          using the Cholesky factorization computed by potrf()
00562     namespace detail {
00563 
00564       template <typename SymmA, typename MatrB>
00565       inline
00566       int potrs (CBLAS_UPLO const uplo, SymmA const& a, MatrB& b) {
00567 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00568         BOOST_STATIC_ASSERT((boost::is_same<
00569           typename traits::matrix_traits<SymmA>::ordering_type,
00570           typename traits::matrix_traits<MatrB>::ordering_type
00571         >::value)); 
00572 #endif 
00573 
00574         CBLAS_ORDER const stor_ord
00575           = enum_cast<CBLAS_ORDER const>
00576           (storage_order<
00577 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00578             typename traits::matrix_traits<SymmA>::ordering_type
00579 #else
00580             typename SymmA::orientation_category 
00581 #endif 
00582            >::value); 
00583 
00584         int const n = traits::matrix_size1 (a);
00585         int const nrhs = stor_ord == CblasColMajor
00586           ? traits::matrix_size2 (b)
00587           : traits::matrix_size1 (b); 
00588         assert (n == traits::matrix_size2 (a)); 
00589         assert (n == (stor_ord == CblasColMajor
00590                       ? traits::matrix_size1 (b)
00591                       : traits::matrix_size2 (b))); 
00592 
00593 #ifndef BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00594         return potrs (stor_ord, uplo, n, nrhs, 
00595 #ifndef BOOST_NO_FUNCTION_TEMPLATE_ORDERING
00596                       traits::matrix_storage (a), 
00597 #else
00598                       traits::matrix_storage_const (a), 
00599 #endif 
00600                       traits::leading_dimension (a),
00601                       traits::matrix_storage (b),
00602                       traits::leading_dimension (b));
00603 #else // BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00604         int ierr; 
00605         if (stor_ord == CblasColMajor)
00606           ierr = potrs (stor_ord, uplo, n, nrhs, 
00607 #ifndef BOOST_NO_FUNCTION_TEMPLATE_ORDERING
00608                         traits::matrix_storage (a), 
00609 #else
00610                         traits::matrix_storage_const (a), 
00611 #endif 
00612                         traits::leading_dimension (a),
00613                         traits::matrix_storage (b),
00614                         traits::leading_dimension (b));
00615         else // ATLAS bug with CblasRowMajor 
00616           ierr = potrs_bug (stor_ord, uplo, n, nrhs, 
00617 #ifndef BOOST_NO_FUNCTION_TEMPLATE_ORDERING
00618                             traits::matrix_storage (a), 
00619 #else
00620                             traits::matrix_storage_const (a), 
00621 #endif 
00622                             traits::leading_dimension (a),
00623                             traits::matrix_storage (b),
00624                             traits::leading_dimension (b));
00625         return ierr; 
00626 #endif // BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00627       }
00628 
00629     } // detail 
00630 
00631     template <typename SymmA, typename MatrB>
00632     inline
00633     int potrs (CBLAS_UPLO const uplo, SymmA const& a, MatrB& b) {
00634 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
00635       BOOST_STATIC_ASSERT((boost::is_same<
00636         typename traits::matrix_traits<SymmA>::matrix_structure, 
00637         traits::general_t
00638       >::value));
00639       BOOST_STATIC_ASSERT((boost::is_same<
00640         typename traits::matrix_traits<MatrB>::matrix_structure, 
00641         traits::general_t
00642       >::value));
00643 #endif
00644       assert (uplo == CblasUpper || uplo == CblasLower); 
00645       return detail::potrs (uplo, a, b); 
00646     }
00647 
00648     template <typename SymmA, typename MatrB>
00649     inline
00650     int potrs (SymmA const& a, MatrB& b) {
00651 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00652       typedef typename traits::matrix_traits<SymmA>::value_type val_t;
00653       BOOST_STATIC_ASSERT((boost::is_same<
00654         typename traits::matrix_traits<SymmA>::matrix_structure, 
00655         typename traits::detail::symm_herm_t<val_t>::type
00656       >::value)); 
00657       BOOST_STATIC_ASSERT((boost::is_same<
00658         typename traits::matrix_traits<MatrB>::matrix_structure, 
00659         traits::general_t
00660       >::value)); 
00661 #endif 
00662 
00663       CBLAS_UPLO const uplo
00664         = enum_cast<CBLAS_UPLO const>
00665         (uplo_triang<
00666 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00667           typename traits::matrix_traits<SymmA>::uplo_type
00668 #else
00669           typename SymmA::packed_category 
00670 #endif 
00671          >::value); 
00672       
00673       return detail::potrs (uplo, a, b); 
00674     }
00675 
00676     template <typename SymmA, typename MatrB>
00677     inline 
00678     int cholesky_substitute (SymmA const& a, MatrB& b) { return potrs (a, b); }
00679 
00680 
00681 #ifdef BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG
00682     // posv(): 'driver' function 
00683     template <typename SymmA, typename MatrB>
00684     inline
00685     int posv (CBLAS_UPLO const uplo, SymmA& a, MatrB& b) {
00686       int ierr = potrf (uplo, a); 
00687       if (ierr == 0)
00688         ierr = potrs (uplo, a, b);
00689       return ierr; 
00690     }
00691 
00692     template <typename SymmA, typename MatrB>
00693     inline
00694     int posv (SymmA& a, MatrB& b) {
00695       int ierr = potrf (a); 
00696       if (ierr == 0)
00697         ierr = potrs (a, b);
00698       return ierr; 
00699     }
00700 
00701     template <typename SymmA, typename MatrB>
00702     inline
00703     int cholesky_solve (SymmA& a, MatrB& b) {
00704       return posv (a, b); 
00705     }
00706 #endif // BOOST_NUMERIC_BINDINGS_ATLAS_POTRF_BUG 
00707 
00708 
00709     // potri(): computes the inverse of a symmetric or Hermitian positive 
00710     //          definite matrix A using the Cholesky factorization 
00711     //          previously computed by potrf() 
00712     namespace detail {
00713 
00714       template <typename SymmA>
00715       inline
00716       int potri (CBLAS_UPLO const uplo, SymmA& a) {
00717 
00718         CBLAS_ORDER const stor_ord
00719           = enum_cast<CBLAS_ORDER const>
00720           (storage_order<
00721 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00722             typename traits::matrix_traits<SymmA>::ordering_type
00723 #else
00724             typename SymmA::orientation_category 
00725 #endif 
00726            >::value); 
00727 
00728         int const n = traits::matrix_size1 (a);
00729         assert (n == traits::matrix_size2 (a)); 
00730 
00731         return potri (stor_ord, uplo, n, 
00732                       traits::matrix_storage (a), 
00733                       traits::leading_dimension (a));
00734       }
00735 
00736     } // detail 
00737 
00738     template <typename SymmA>
00739     inline
00740     int potri (CBLAS_UPLO const uplo, SymmA& a) {
00741 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK
00742       BOOST_STATIC_ASSERT((boost::is_same<
00743         typename traits::matrix_traits<SymmA>::matrix_structure, 
00744         traits::general_t
00745       >::value));
00746 #endif
00747       assert (uplo == CblasUpper || uplo == CblasLower); 
00748       return detail::potri (uplo, a); 
00749     } 
00750 
00751     template <typename SymmA>
00752     inline
00753     int potri (SymmA& a) {
00754 #ifndef BOOST_NUMERIC_BINDINGS_NO_STRUCTURE_CHECK 
00755       typedef typename traits::matrix_traits<SymmA>::value_type val_t;
00756       BOOST_STATIC_ASSERT((boost::is_same<
00757         typename traits::matrix_traits<SymmA>::matrix_structure, 
00758         typename traits::detail::symm_herm_t<val_t>::type
00759       >::value)); 
00760 #endif 
00761 
00762       CBLAS_UPLO const uplo
00763         = enum_cast<CBLAS_UPLO const>
00764         (uplo_triang<
00765 #ifndef BOOST_NUMERIC_BINDINGS_POOR_MANS_TRAITS
00766           typename traits::matrix_traits<SymmA>::uplo_type
00767 #else
00768           typename SymmA::packed_category 
00769 #endif 
00770          >::value); 
00771       
00772       return detail::potri (uplo, a); 
00773     }
00774 
00775     template <typename SymmA>
00776     inline
00777     int cholesky_invert (SymmA& a) { return potri (a); }
00778 
00779 
00780   } // namespace atlas
00781 
00782 }}} 
00783 
00784 #endif // BOOST_NUMERIC_BINDINGS_CLAPACK_HPP

Generated on Wed Nov 23 19:00:00 2011 for FreeCAD by  doxygen 1.6.1