123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476 |
- // Copyright 2008-2016 Conrad Sanderson (http://conradsanderson.id.au)
- // Copyright 2008-2016 National ICT Australia (NICTA)
- //
- // Licensed under the Apache License, Version 2.0 (the "License");
- // you may not use this file except in compliance with the License.
- // You may obtain a copy of the License at
- // http://www.apache.org/licenses/LICENSE-2.0
- //
- // Unless required by applicable law or agreed to in writing, software
- // distributed under the License is distributed on an "AS IS" BASIS,
- // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
- // See the License for the specific language governing permissions and
- // limitations under the License.
- // ------------------------------------------------------------------------
- //! \addtogroup auxlib
- //! @{
- //! interface functions for accessing decompositions in LAPACK and ATLAS
- class auxlib
- {
- public:
-
-
- template<const uword row, const uword col>
- struct pos
- {
- static const uword n2 = row + col*2;
- static const uword n3 = row + col*3;
- static const uword n4 = row + col*4;
- };
-
-
- //
- // inv
-
- template<typename eT>
- inline static bool inv(Mat<eT>& out, const Mat<eT>& A);
-
- template<typename eT>
- arma_cold inline static bool inv_tiny(Mat<eT>& out, const Mat<eT>& X);
-
- template<typename eT, typename T1>
- inline static bool inv_tr(Mat<eT>& out, const Base<eT,T1>& X, const uword layout);
-
- template<typename eT, typename T1>
- inline static bool inv_sympd(Mat<eT>& out, const Base<eT,T1>& X);
-
- template<typename eT>
- arma_cold inline static bool inv_sympd_tiny(Mat<eT>& out, const Mat<eT>& X);
-
-
-
- //
- // det
-
- template<typename eT>
- inline static eT det(const Mat<eT>& A);
-
- template<typename eT>
- arma_cold inline static eT det_tinymat(const Mat<eT>& X, const uword N);
-
- template<typename eT>
- inline static eT det_lapack(const Mat<eT>& X);
-
-
- //
- // log_det
-
- template<typename eT, typename T1>
- inline static bool log_det(eT& out_val, typename get_pod_type<eT>::result& out_sign, const Base<eT,T1>& X);
-
-
- //
- // lu
-
- template<typename eT, typename T1>
- inline static bool lu(Mat<eT>& L, Mat<eT>& U, podarray<blas_int>& ipiv, const Base<eT,T1>& X);
-
- template<typename eT, typename T1>
- inline static bool lu(Mat<eT>& L, Mat<eT>& U, Mat<eT>& P, const Base<eT,T1>& X);
-
- template<typename eT, typename T1>
- inline static bool lu(Mat<eT>& L, Mat<eT>& U, const Base<eT,T1>& X);
-
-
- //
- // eig_gen
-
- template<typename T1>
- inline static bool eig_gen(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& vecs, const bool vecs_on, const Base<typename T1::pod_type,T1>& expr);
-
- template<typename T1>
- inline static bool eig_gen(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& vecs, const bool vecs_on, const Base< std::complex<typename T1::pod_type>, T1 >& expr);
-
-
- //
- // eig_gen_balance
-
- template<typename T1>
- inline static bool eig_gen_balance(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& vecs, const bool vecs_on, const Base<typename T1::pod_type,T1>& expr);
-
- template<typename T1>
- inline static bool eig_gen_balance(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& vecs, const bool vecs_on, const Base< std::complex<typename T1::pod_type>, T1 >& expr);
-
-
- //
- // eig_gen_twosided
-
- template<typename T1>
- inline static bool eig_gen_twosided(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& lvecs, Mat< std::complex<typename T1::pod_type> >& rvecs, const Base<typename T1::pod_type,T1>& expr);
-
- template<typename T1>
- inline static bool eig_gen_twosided(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& lvecs, Mat< std::complex<typename T1::pod_type> >& rvecs, const Base< std::complex<typename T1::pod_type>, T1 >& expr);
-
-
- //
- // eig_gen_twosided_balance
-
- template<typename T1>
- inline static bool eig_gen_twosided_balance(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& lvecs, Mat< std::complex<typename T1::pod_type> >& rvecs, const Base<typename T1::pod_type,T1>& expr);
-
- template<typename T1>
- inline static bool eig_gen_twosided_balance(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& lvecs, Mat< std::complex<typename T1::pod_type> >& rvecs, const Base< std::complex<typename T1::pod_type>, T1 >& expr);
-
-
- //
- // eig_pair
-
- template<typename T1, typename T2>
- inline static bool eig_pair(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& vecs, const bool vecs_on, const Base<typename T1::pod_type,T1>& A_expr, const Base<typename T1::pod_type,T2>& B_expr);
-
- template<typename T1, typename T2>
- inline static bool eig_pair(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& vecs, const bool vecs_on, const Base< std::complex<typename T1::pod_type>, T1 >& A_expr, const Base< std::complex<typename T1::pod_type>, T2 >& B_expr);
-
-
- //
- // eig_pair_twosided
-
- template<typename T1, typename T2>
- inline static bool eig_pair_twosided(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& lvecs, Mat< std::complex<typename T1::pod_type> >& rvecs, const Base<typename T1::pod_type,T1>& A_expr, const Base<typename T1::pod_type,T2>& B_expr);
-
- template<typename T1, typename T2>
- inline static bool eig_pair_twosided(Mat< std::complex<typename T1::pod_type> >& vals, Mat< std::complex<typename T1::pod_type> >& lvecs, Mat< std::complex<typename T1::pod_type> >& rvecs, const Base< std::complex<typename T1::pod_type>, T1 >& A_expr, const Base< std::complex<typename T1::pod_type>, T2 >& B_expr);
-
-
- //
- // eig_sym
-
- template<typename eT, typename T1>
- inline static bool eig_sym(Col<eT>& eigval, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool eig_sym(Col<T>& eigval, const Base<std::complex<T>,T1>& X);
-
- template<typename eT>
- inline static bool eig_sym(Col<eT>& eigval, Mat<eT>& eigvec, const Mat<eT>& X);
-
- template<typename T>
- inline static bool eig_sym(Col<T>& eigval, Mat< std::complex<T> >& eigvec, const Mat< std::complex<T> >& X);
-
- template<typename eT>
- inline static bool eig_sym_dc(Col<eT>& eigval, Mat<eT>& eigvec, const Mat<eT>& X);
-
- template<typename T>
- inline static bool eig_sym_dc(Col<T>& eigval, Mat< std::complex<T> >& eigvec, const Mat< std::complex<T> >& X);
-
-
- //
- // chol
-
- template<typename eT>
- inline static bool chol_simple(Mat<eT>& X);
-
- template<typename eT>
- inline static bool chol(Mat<eT>& X, const uword layout);
-
- template<typename eT>
- inline static bool chol_band(Mat<eT>& X, const uword KD, const uword layout);
-
- template<typename T>
- inline static bool chol_band(Mat< std::complex<T> >& X, const uword KD, const uword layout);
-
- template<typename eT>
- inline static bool chol_band_common(Mat<eT>& X, const uword KD, const uword layout);
-
-
- //
- // hessenberg decomposition
-
- template<typename eT, typename T1>
- inline static bool hess(Mat<eT>& H, const Base<eT,T1>& X, Col<eT>& tao);
-
-
- //
- // qr
-
- template<typename eT, typename T1>
- inline static bool qr(Mat<eT>& Q, Mat<eT>& R, const Base<eT,T1>& X);
-
- template<typename eT, typename T1>
- inline static bool qr_econ(Mat<eT>& Q, Mat<eT>& R, const Base<eT,T1>& X);
-
- template<typename eT, typename T1>
- inline static bool qr_pivot(Mat<eT>& Q, Mat<eT>& R, Mat<uword>& P, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool qr_pivot(Mat< std::complex<T> >& Q, Mat< std::complex<T> >& R, Mat<uword>& P, const Base<std::complex<T>,T1>& X);
-
-
- //
- // svd
-
- template<typename eT, typename T1>
- inline static bool svd(Col<eT>& S, const Base<eT,T1>& X, uword& n_rows, uword& n_cols);
-
- template<typename T, typename T1>
- inline static bool svd(Col<T>& S, const Base<std::complex<T>, T1>& X, uword& n_rows, uword& n_cols);
-
- template<typename eT, typename T1>
- inline static bool svd(Col<eT>& S, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool svd(Col<T>& S, const Base<std::complex<T>, T1>& X);
-
- template<typename eT, typename T1>
- inline static bool svd(Mat<eT>& U, Col<eT>& S, Mat<eT>& V, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool svd(Mat< std::complex<T> >& U, Col<T>& S, Mat< std::complex<T> >& V, const Base< std::complex<T>, T1>& X);
-
- template<typename eT, typename T1>
- inline static bool svd_econ(Mat<eT>& U, Col<eT>& S, Mat<eT>& V, const Base<eT,T1>& X, const char mode);
-
- template<typename T, typename T1>
- inline static bool svd_econ(Mat< std::complex<T> >& U, Col<T>& S, Mat< std::complex<T> >& V, const Base< std::complex<T>, T1>& X, const char mode);
-
-
- template<typename eT, typename T1>
- inline static bool svd_dc(Col<eT>& S, const Base<eT,T1>& X, uword& n_rows, uword& n_cols);
-
- template<typename T, typename T1>
- inline static bool svd_dc(Col<T>& S, const Base<std::complex<T>, T1>& X, uword& n_rows, uword& n_cols);
-
- template<typename eT, typename T1>
- inline static bool svd_dc(Col<eT>& S, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool svd_dc(Col<T>& S, const Base<std::complex<T>, T1>& X);
-
-
- template<typename eT, typename T1>
- inline static bool svd_dc(Mat<eT>& U, Col<eT>& S, Mat<eT>& V, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool svd_dc(Mat< std::complex<T> >& U, Col<T>& S, Mat< std::complex<T> >& V, const Base< std::complex<T>, T1>& X);
-
- template<typename eT, typename T1>
- inline static bool svd_dc_econ(Mat<eT>& U, Col<eT>& S, Mat<eT>& V, const Base<eT,T1>& X);
-
- template<typename T, typename T1>
- inline static bool svd_dc_econ(Mat< std::complex<T> >& U, Col<T>& S, Mat< std::complex<T> >& V, const Base< std::complex<T>, T1>& X);
-
-
- //
- // solve
-
- template<typename T1>
- arma_cold inline static bool solve_square_tiny(Mat<typename T1::elem_type>& out, const Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_square_fast(Mat<typename T1::elem_type>& out, Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_square_rcond(Mat<typename T1::elem_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_square_refine(Mat<typename T1::pod_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::pod_type>& A, const Base<typename T1::pod_type,T1>& B_expr, const bool equilibrate, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_square_refine(Mat< std::complex<typename T1::pod_type> >& out, typename T1::pod_type& out_rcond, Mat< std::complex<typename T1::pod_type> >& A, const Base<std::complex<typename T1::pod_type>,T1>& B_expr, const bool equilibrate, const bool allow_ugly);
-
- //
-
- template<typename T1>
- inline static bool solve_sympd_fast(Mat<typename T1::elem_type>& out, Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_sympd_fast_common(Mat<typename T1::elem_type>& out, Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_sympd_rcond(Mat<typename T1::pod_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::pod_type>& A, const Base<typename T1::pod_type,T1>& B_expr, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_sympd_rcond(Mat< std::complex<typename T1::pod_type> >& out, typename T1::pod_type& out_rcond, Mat< std::complex<typename T1::pod_type> >& A, const Base< std::complex<typename T1::pod_type>,T1>& B_expr, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_sympd_refine(Mat<typename T1::pod_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::pod_type>& A, const Base<typename T1::pod_type,T1>& B_expr, const bool equilibrate, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_sympd_refine(Mat< std::complex<typename T1::pod_type> >& out, typename T1::pod_type& out_rcond, Mat< std::complex<typename T1::pod_type> >& A, const Base<std::complex<typename T1::pod_type>,T1>& B_expr, const bool equilibrate, const bool allow_ugly);
-
- //
-
- template<typename T1>
- inline static bool solve_rect_fast(Mat<typename T1::elem_type>& out, Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_rect_rcond(Mat<typename T1::elem_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr, const bool allow_ugly);
-
- //
-
- template<typename T1>
- inline static bool solve_approx_svd(Mat<typename T1::pod_type>& out, Mat<typename T1::pod_type>& A, const Base<typename T1::pod_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_approx_svd(Mat< std::complex<typename T1::pod_type> >& out, Mat< std::complex<typename T1::pod_type> >& A, const Base<std::complex<typename T1::pod_type>,T1>& B_expr);
-
- //
-
- template<typename T1>
- inline static bool solve_trimat_fast(Mat<typename T1::elem_type>& out, const Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr, const uword layout);
-
- template<typename T1>
- inline static bool solve_trimat_rcond(Mat<typename T1::elem_type>& out, typename T1::pod_type& out_rcond, const Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr, const uword layout, const bool allow_ugly);
-
- //
-
- template<typename T1>
- inline static bool solve_band_fast(Mat<typename T1::pod_type>& out, Mat<typename T1::pod_type>& A, const uword KL, const uword KU, const Base<typename T1::pod_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_band_fast(Mat< std::complex<typename T1::pod_type> >& out, Mat< std::complex<typename T1::pod_type> >& A, const uword KL, const uword KU, const Base< std::complex<typename T1::pod_type>,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_band_fast_common(Mat<typename T1::elem_type>& out, const Mat<typename T1::elem_type>& A, const uword KL, const uword KU, const Base<typename T1::elem_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_band_rcond(Mat<typename T1::pod_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::pod_type>& A, const uword KL, const uword KU, const Base<typename T1::pod_type,T1>& B_expr, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_band_rcond(Mat< std::complex<typename T1::pod_type> >& out, typename T1::pod_type& out_rcond, Mat< std::complex<typename T1::pod_type> >& A, const uword KL, const uword KU, const Base< std::complex<typename T1::pod_type>,T1>& B_expr, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_band_rcond_common(Mat<typename T1::elem_type>& out, typename T1::pod_type& out_rcond, const Mat<typename T1::elem_type>& A, const uword KL, const uword KU, const Base<typename T1::elem_type,T1>& B_expr, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_band_refine(Mat<typename T1::pod_type>& out, typename T1::pod_type& out_rcond, Mat<typename T1::pod_type>& A, const uword KL, const uword KU, const Base<typename T1::pod_type,T1>& B_expr, const bool equilibrate, const bool allow_ugly);
-
- template<typename T1>
- inline static bool solve_band_refine(Mat< std::complex<typename T1::pod_type> >& out, typename T1::pod_type& out_rcond, Mat< std::complex<typename T1::pod_type> >& A, const uword KL, const uword KU, const Base<std::complex<typename T1::pod_type>,T1>& B_expr, const bool equilibrate, const bool allow_ugly);
-
- //
-
- template<typename T1>
- inline static bool solve_tridiag_fast(Mat<typename T1::pod_type>& out, Mat<typename T1::pod_type>& A, const Base<typename T1::pod_type,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_tridiag_fast(Mat< std::complex<typename T1::pod_type> >& out, Mat< std::complex<typename T1::pod_type> >& A, const Base< std::complex<typename T1::pod_type>,T1>& B_expr);
-
- template<typename T1>
- inline static bool solve_tridiag_fast_common(Mat<typename T1::elem_type>& out, const Mat<typename T1::elem_type>& A, const Base<typename T1::elem_type,T1>& B_expr);
-
-
- //
- // Schur decomposition
-
- template<typename eT, typename T1>
- inline static bool schur(Mat<eT>& U, Mat<eT>& S, const Base<eT,T1>& X, const bool calc_U = true);
-
- template<typename T, typename T1>
- inline static bool schur(Mat<std::complex<T> >& U, Mat<std::complex<T> >& S, const Base<std::complex<T>,T1>& X, const bool calc_U = true);
-
- template<typename T>
- inline static bool schur(Mat<std::complex<T> >& U, Mat<std::complex<T> >& S, const bool calc_U = true);
-
- //
- // syl (solution of the Sylvester equation AX + XB = C)
-
- template<typename eT>
- inline static bool syl(Mat<eT>& X, const Mat<eT>& A, const Mat<eT>& B, const Mat<eT>& C);
-
-
- //
- // QZ decomposition
-
- template<typename T, typename T1, typename T2>
- inline static bool qz(Mat<T>& A, Mat<T>& B, Mat<T>& vsl, Mat<T>& vsr, const Base<T,T1>& X_expr, const Base<T,T2>& Y_expr, const char mode);
-
- template<typename T, typename T1, typename T2>
- inline static bool qz(Mat< std::complex<T> >& A, Mat< std::complex<T> >& B, Mat< std::complex<T> >& vsl, Mat< std::complex<T> >& vsr, const Base< std::complex<T>, T1 >& X_expr, const Base< std::complex<T>, T2 >& Y_expr, const char mode);
-
-
- //
- // rcond
-
- template<typename eT>
- inline static eT rcond(Mat<eT>& A);
-
- template<typename T>
- inline static T rcond(Mat< std::complex<T> >& A);
-
- template<typename eT>
- inline static eT rcond_sympd(Mat<eT>& A, bool& calc_ok);
-
- template<typename T>
- inline static T rcond_sympd(Mat< std::complex<T> >& A, bool& calc_ok);
-
- template<typename eT>
- inline static eT rcond_trimat(const Mat<eT>& A, const uword layout);
-
- template<typename T>
- inline static T rcond_trimat(const Mat< std::complex<T> >& A, const uword layout);
-
-
- //
- // lu_rcond (rcond from pre-computed LU decomposition)
-
- template<typename eT>
- inline static eT lu_rcond(const Mat<eT>& A, const eT norm_val);
-
- template<typename T>
- inline static T lu_rcond(const Mat< std::complex<T> >& A, const T norm_val);
-
- template<typename eT>
- inline static eT lu_rcond_sympd(const Mat<eT>& A, const eT norm_val);
-
- template<typename T>
- inline static T lu_rcond_sympd(const Mat< std::complex<T> >& A, const T norm_val);
-
- template<typename eT>
- inline static eT lu_rcond_band(const Mat<eT>& AB, const uword KL, const uword KU, const podarray<blas_int>& ipiv, const eT norm_val);
-
- template<typename T>
- inline static T lu_rcond_band(const Mat< std::complex<T> >& AB, const uword KL, const uword KU, const podarray<blas_int>& ipiv, const T norm_val);
-
-
- //
- // misc
-
- template<typename T1>
- inline static bool crippled_lapack(const Base<typename T1::elem_type, T1>&);
-
- template<typename T1>
- inline static typename T1::pod_type epsilon_lapack(const Base<typename T1::elem_type, T1>&);
-
- template<typename eT>
- inline static bool rudimentary_sym_check(const Mat<eT>& X);
-
- template<typename T>
- inline static bool rudimentary_sym_check(const Mat< std::complex<T> >& X);
- };
- namespace qz_helper
- {
- template<typename T> inline blas_int select_lhp(const T* x_ptr, const T* y_ptr, const T* z_ptr);
- template<typename T> inline blas_int select_rhp(const T* x_ptr, const T* y_ptr, const T* z_ptr);
- template<typename T> inline blas_int select_iuc(const T* x_ptr, const T* y_ptr, const T* z_ptr);
- template<typename T> inline blas_int select_ouc(const T* x_ptr, const T* y_ptr, const T* z_ptr);
-
- template<typename T> inline blas_int cx_select_lhp(const std::complex<T>* x_ptr, const std::complex<T>* y_ptr);
- template<typename T> inline blas_int cx_select_rhp(const std::complex<T>* x_ptr, const std::complex<T>* y_ptr);
- template<typename T> inline blas_int cx_select_iuc(const std::complex<T>* x_ptr, const std::complex<T>* y_ptr);
- template<typename T> inline blas_int cx_select_ouc(const std::complex<T>* x_ptr, const std::complex<T>* y_ptr);
-
- template<typename T> inline void_ptr ptr_cast(blas_int (*function)(const T*, const T*, const T*));
- template<typename T> inline void_ptr ptr_cast(blas_int (*function)(const std::complex<T>*, const std::complex<T>*));
- }
- //! @}
|