Mercurial > octave-nkf
view liboctave/numeric/CmplxSVD.cc @ 17769:49a5a4be04a1
maint: Use GNU style coding conventions for code in liboctave/
* liboctave/array/Array-C.cc, liboctave/array/Array-b.cc,
liboctave/array/Array-ch.cc, liboctave/array/Array-d.cc,
liboctave/array/Array-f.cc, liboctave/array/Array-fC.cc,
liboctave/array/Array-util.cc, liboctave/array/Array-util.h,
liboctave/array/Array.cc, liboctave/array/Array.h, liboctave/array/Array3.h,
liboctave/array/CColVector.cc, liboctave/array/CColVector.h,
liboctave/array/CDiagMatrix.cc, liboctave/array/CDiagMatrix.h,
liboctave/array/CMatrix.cc, liboctave/array/CMatrix.h,
liboctave/array/CNDArray.cc, liboctave/array/CNDArray.h,
liboctave/array/CRowVector.cc, liboctave/array/CRowVector.h,
liboctave/array/CSparse.cc, liboctave/array/CSparse.h,
liboctave/array/DiagArray2.h, liboctave/array/MArray.cc,
liboctave/array/MArray.h, liboctave/array/MDiagArray2.cc,
liboctave/array/MDiagArray2.h, liboctave/array/MSparse.cc,
liboctave/array/MSparse.h, liboctave/array/MatrixType.cc,
liboctave/array/MatrixType.h, liboctave/array/PermMatrix.h,
liboctave/array/Range.cc, liboctave/array/Range.h, liboctave/array/Sparse.cc,
liboctave/array/Sparse.h, liboctave/array/boolMatrix.cc,
liboctave/array/boolMatrix.h, liboctave/array/boolNDArray.cc,
liboctave/array/boolNDArray.h, liboctave/array/boolSparse.cc,
liboctave/array/boolSparse.h, liboctave/array/chMatrix.cc,
liboctave/array/chMatrix.h, liboctave/array/chNDArray.cc,
liboctave/array/chNDArray.h, liboctave/array/dColVector.h,
liboctave/array/dDiagMatrix.cc, liboctave/array/dDiagMatrix.h,
liboctave/array/dMatrix.cc, liboctave/array/dMatrix.h,
liboctave/array/dNDArray.cc, liboctave/array/dNDArray.h,
liboctave/array/dRowVector.h, liboctave/array/dSparse.cc,
liboctave/array/dSparse.h, liboctave/array/dim-vector.cc,
liboctave/array/dim-vector.h, liboctave/array/fCColVector.cc,
liboctave/array/fCColVector.h, liboctave/array/fCDiagMatrix.cc,
liboctave/array/fCDiagMatrix.h, liboctave/array/fCMatrix.cc,
liboctave/array/fCMatrix.h, liboctave/array/fCNDArray.cc,
liboctave/array/fCNDArray.h, liboctave/array/fCRowVector.cc,
liboctave/array/fCRowVector.h, liboctave/array/fColVector.h,
liboctave/array/fDiagMatrix.cc, liboctave/array/fDiagMatrix.h,
liboctave/array/fMatrix.cc, liboctave/array/fMatrix.h,
liboctave/array/fNDArray.cc, liboctave/array/fNDArray.h,
liboctave/array/fRowVector.h, liboctave/array/idx-vector.cc,
liboctave/array/idx-vector.h, liboctave/array/intNDArray.cc,
liboctave/array/intNDArray.h, liboctave/cruft/misc/blaswrap.c,
liboctave/cruft/misc/quit.cc, liboctave/numeric/CmplxCHOL.cc,
liboctave/numeric/CmplxCHOL.h, liboctave/numeric/CmplxGEPBAL.cc,
liboctave/numeric/CmplxGEPBAL.h, liboctave/numeric/CmplxHESS.h,
liboctave/numeric/CmplxLU.cc, liboctave/numeric/CmplxLU.h,
liboctave/numeric/CmplxQR.cc, liboctave/numeric/CmplxQRP.cc,
liboctave/numeric/CmplxQRP.h, liboctave/numeric/CmplxSCHUR.h,
liboctave/numeric/CmplxSVD.cc, liboctave/numeric/CmplxSVD.h,
liboctave/numeric/CollocWt.h, liboctave/numeric/DAE.h,
liboctave/numeric/DAEFunc.h, liboctave/numeric/DAERT.h,
liboctave/numeric/DAERTFunc.h, liboctave/numeric/DASPK.cc,
liboctave/numeric/DASRT.cc, liboctave/numeric/DASRT.h,
liboctave/numeric/DASSL.cc, liboctave/numeric/DET.h, liboctave/numeric/EIG.cc,
liboctave/numeric/EIG.h, liboctave/numeric/LSODE.cc, liboctave/numeric/ODE.h,
liboctave/numeric/ODEFunc.h, liboctave/numeric/ODES.h,
liboctave/numeric/ODESFunc.h, liboctave/numeric/Quad.cc,
liboctave/numeric/Quad.h, liboctave/numeric/SparseCmplxCHOL.h,
liboctave/numeric/SparseCmplxLU.cc, liboctave/numeric/SparseCmplxLU.h,
liboctave/numeric/SparseCmplxQR.cc, liboctave/numeric/SparseCmplxQR.h,
liboctave/numeric/SparseQR.cc, liboctave/numeric/SparseQR.h,
liboctave/numeric/SparsedbleCHOL.h, liboctave/numeric/SparsedbleLU.cc,
liboctave/numeric/SparsedbleLU.h, liboctave/numeric/base-aepbal.h,
liboctave/numeric/base-dae.h, liboctave/numeric/base-de.h,
liboctave/numeric/base-lu.cc, liboctave/numeric/base-lu.h,
liboctave/numeric/base-min.h, liboctave/numeric/base-qr.h,
liboctave/numeric/bsxfun.h, liboctave/numeric/dbleCHOL.cc,
liboctave/numeric/dbleCHOL.h, liboctave/numeric/dbleGEPBAL.h,
liboctave/numeric/dbleHESS.h, liboctave/numeric/dbleLU.cc,
liboctave/numeric/dbleLU.h, liboctave/numeric/dbleQR.cc,
liboctave/numeric/dbleQRP.cc, liboctave/numeric/dbleQRP.h,
liboctave/numeric/dbleSCHUR.cc, liboctave/numeric/dbleSCHUR.h,
liboctave/numeric/dbleSVD.cc, liboctave/numeric/dbleSVD.h,
liboctave/numeric/eigs-base.cc, liboctave/numeric/fCmplxAEPBAL.cc,
liboctave/numeric/fCmplxAEPBAL.h, liboctave/numeric/fCmplxCHOL.cc,
liboctave/numeric/fCmplxCHOL.h, liboctave/numeric/fCmplxGEPBAL.cc,
liboctave/numeric/fCmplxGEPBAL.h, liboctave/numeric/fCmplxHESS.h,
liboctave/numeric/fCmplxLU.cc, liboctave/numeric/fCmplxLU.h,
liboctave/numeric/fCmplxQR.cc, liboctave/numeric/fCmplxQR.h,
liboctave/numeric/fCmplxQRP.cc, liboctave/numeric/fCmplxQRP.h,
liboctave/numeric/fCmplxSCHUR.cc, liboctave/numeric/fCmplxSCHUR.h,
liboctave/numeric/fCmplxSVD.h, liboctave/numeric/fEIG.cc,
liboctave/numeric/fEIG.h, liboctave/numeric/floatCHOL.cc,
liboctave/numeric/floatCHOL.h, liboctave/numeric/floatGEPBAL.cc,
liboctave/numeric/floatGEPBAL.h, liboctave/numeric/floatHESS.h,
liboctave/numeric/floatLU.cc, liboctave/numeric/floatLU.h,
liboctave/numeric/floatQR.cc, liboctave/numeric/floatQRP.cc,
liboctave/numeric/floatQRP.h, liboctave/numeric/floatSCHUR.cc,
liboctave/numeric/floatSCHUR.h, liboctave/numeric/floatSVD.cc,
liboctave/numeric/floatSVD.h, liboctave/numeric/lo-mappers.cc,
liboctave/numeric/lo-mappers.h, liboctave/numeric/lo-specfun.cc,
liboctave/numeric/lo-specfun.h, liboctave/numeric/oct-convn.cc,
liboctave/numeric/oct-fftw.cc, liboctave/numeric/oct-fftw.h,
liboctave/numeric/oct-norm.cc, liboctave/numeric/oct-rand.cc,
liboctave/numeric/oct-rand.h, liboctave/numeric/randgamma.c,
liboctave/numeric/randgamma.h, liboctave/numeric/randmtzig.c,
liboctave/numeric/randpoisson.c, liboctave/numeric/randpoisson.h,
liboctave/numeric/sparse-base-chol.h, liboctave/numeric/sparse-base-lu.h,
liboctave/numeric/sparse-dmsolve.cc, liboctave/operators/Sparse-diag-op-defs.h,
liboctave/operators/Sparse-op-defs.h, liboctave/operators/mx-inlines.cc,
liboctave/system/dir-ops.h, liboctave/system/file-ops.cc,
liboctave/system/file-stat.cc, liboctave/system/file-stat.h,
liboctave/system/lo-sysdep.cc, liboctave/system/lo-sysdep.h,
liboctave/system/mach-info.cc, liboctave/system/mach-info.h,
liboctave/system/oct-env.cc, liboctave/system/oct-group.cc,
liboctave/system/oct-syscalls.cc, liboctave/system/oct-syscalls.h,
liboctave/system/oct-time.h, liboctave/system/tempname.c,
liboctave/util/action-container.h, liboctave/util/base-list.h,
liboctave/util/cmd-edit.cc, liboctave/util/cmd-edit.h,
liboctave/util/cmd-hist.cc, liboctave/util/cmd-hist.h,
liboctave/util/data-conv.cc, liboctave/util/data-conv.h,
liboctave/util/kpse.cc, liboctave/util/lo-array-gripes.cc,
liboctave/util/lo-cieee.c, liboctave/util/lo-regexp.cc,
liboctave/util/lo-utils.cc, liboctave/util/oct-alloc.cc,
liboctave/util/oct-base64.cc, liboctave/util/oct-binmap.h,
liboctave/util/oct-cmplx.h, liboctave/util/oct-glob.cc,
liboctave/util/oct-inttypes.cc, liboctave/util/oct-inttypes.h,
liboctave/util/oct-locbuf.cc, liboctave/util/oct-locbuf.h,
liboctave/util/oct-mem.h, liboctave/util/oct-mutex.cc,
liboctave/util/oct-refcount.h, liboctave/util/oct-shlib.cc,
liboctave/util/oct-shlib.h, liboctave/util/oct-sort.cc,
liboctave/util/oct-sort.h, liboctave/util/pathsearch.cc,
liboctave/util/pathsearch.h, liboctave/util/sparse-util.cc,
liboctave/util/str-vec.cc, liboctave/util/str-vec.h,
liboctave/util/unwind-prot.h, liboctave/util/url-transfer.cc,
liboctave/util/url-transfer.h: Use GNU style coding conventions.
author | Rik <rik@octave.org> |
---|---|
date | Sat, 26 Oct 2013 18:57:05 -0700 |
parents | d63878346099 |
children | 4197fc428c7d |
line wrap: on
line source
/* Copyright (C) 1994-2013 John W. Eaton This file is part of Octave. Octave is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version. Octave is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with Octave; see the file COPYING. If not, see <http://www.gnu.org/licenses/>. */ #ifdef HAVE_CONFIG_H #include <config.h> #endif #include "CmplxSVD.h" #include "f77-fcn.h" #include "lo-error.h" #include "oct-locbuf.h" extern "C" { F77_RET_T F77_FUNC (zgesvd, ZGESVD) (F77_CONST_CHAR_ARG_DECL, F77_CONST_CHAR_ARG_DECL, const octave_idx_type&, const octave_idx_type&, Complex*, const octave_idx_type&, double*, Complex*, const octave_idx_type&, Complex*, const octave_idx_type&, Complex*, const octave_idx_type&, double*, octave_idx_type& F77_CHAR_ARG_LEN_DECL F77_CHAR_ARG_LEN_DECL); F77_RET_T F77_FUNC (zgesdd, ZGESDD) (F77_CONST_CHAR_ARG_DECL, const octave_idx_type&, const octave_idx_type&, Complex*, const octave_idx_type&, double*, Complex*, const octave_idx_type&, Complex*, const octave_idx_type&, Complex*, const octave_idx_type&, double*, octave_idx_type *, octave_idx_type& F77_CHAR_ARG_LEN_DECL); } ComplexMatrix ComplexSVD::left_singular_matrix (void) const { if (type_computed == SVD::sigma_only) { (*current_liboctave_error_handler) ("ComplexSVD: U not computed because type == SVD::sigma_only"); return ComplexMatrix (); } else return left_sm; } ComplexMatrix ComplexSVD::right_singular_matrix (void) const { if (type_computed == SVD::sigma_only) { (*current_liboctave_error_handler) ("ComplexSVD: V not computed because type == SVD::sigma_only"); return ComplexMatrix (); } else return right_sm; } octave_idx_type ComplexSVD::init (const ComplexMatrix& a, SVD::type svd_type, SVD::driver svd_driver) { octave_idx_type info; octave_idx_type m = a.rows (); octave_idx_type n = a.cols (); ComplexMatrix atmp = a; Complex *tmp_data = atmp.fortran_vec (); octave_idx_type min_mn = m < n ? m : n; octave_idx_type max_mn = m > n ? m : n; char jobu = 'A'; char jobv = 'A'; octave_idx_type ncol_u = m; octave_idx_type nrow_vt = n; octave_idx_type nrow_s = m; octave_idx_type ncol_s = n; switch (svd_type) { case SVD::economy: jobu = jobv = 'S'; ncol_u = nrow_vt = nrow_s = ncol_s = min_mn; break; case SVD::sigma_only: // Note: for this case, both jobu and jobv should be 'N', but // there seems to be a bug in dgesvd from Lapack V2.0. To // demonstrate the bug, set both jobu and jobv to 'N' and find // the singular values of [eye(3), eye(3)]. The result is // [-sqrt(2), -sqrt(2), -sqrt(2)]. // // For Lapack 3.0, this problem seems to be fixed. jobu = jobv = 'N'; ncol_u = nrow_vt = 1; break; default: break; } type_computed = svd_type; if (! (jobu == 'N' || jobu == 'O')) left_sm.resize (m, ncol_u); Complex *u = left_sm.fortran_vec (); sigma.resize (nrow_s, ncol_s); double *s_vec = sigma.fortran_vec (); if (! (jobv == 'N' || jobv == 'O')) right_sm.resize (nrow_vt, n); Complex *vt = right_sm.fortran_vec (); // Query ZGESVD for the correct dimension of WORK. octave_idx_type lwork = -1; Array<Complex> work (dim_vector (1, 1)); octave_idx_type one = 1; octave_idx_type m1 = std::max (m, one); octave_idx_type nrow_vt1 = std::max (nrow_vt, one); if (svd_driver == SVD::GESVD) { octave_idx_type lrwork = 5*max_mn; Array<double> rwork (dim_vector (lrwork, 1)); F77_XFCN (zgesvd, ZGESVD, (F77_CONST_CHAR_ARG2 (&jobu, 1), F77_CONST_CHAR_ARG2 (&jobv, 1), m, n, tmp_data, m1, s_vec, u, m1, vt, nrow_vt1, work.fortran_vec (), lwork, rwork.fortran_vec (), info F77_CHAR_ARG_LEN (1) F77_CHAR_ARG_LEN (1))); lwork = static_cast<octave_idx_type> (work(0).real ()); work.resize (dim_vector (lwork, 1)); F77_XFCN (zgesvd, ZGESVD, (F77_CONST_CHAR_ARG2 (&jobu, 1), F77_CONST_CHAR_ARG2 (&jobv, 1), m, n, tmp_data, m1, s_vec, u, m1, vt, nrow_vt1, work.fortran_vec (), lwork, rwork.fortran_vec (), info F77_CHAR_ARG_LEN (1) F77_CHAR_ARG_LEN (1))); } else if (svd_driver == SVD::GESDD) { assert (jobu == jobv); char jobz = jobu; octave_idx_type lrwork; if (jobz == 'N') lrwork = 7*min_mn; else lrwork = 5*min_mn*min_mn + 5*min_mn; Array<double> rwork (dim_vector (lrwork, 1)); OCTAVE_LOCAL_BUFFER (octave_idx_type, iwork, 8*min_mn); F77_XFCN (zgesdd, ZGESDD, (F77_CONST_CHAR_ARG2 (&jobz, 1), m, n, tmp_data, m1, s_vec, u, m1, vt, nrow_vt1, work.fortran_vec (), lwork, rwork.fortran_vec (), iwork, info F77_CHAR_ARG_LEN (1))); lwork = static_cast<octave_idx_type> (work(0).real ()); work.resize (dim_vector (lwork, 1)); F77_XFCN (zgesdd, ZGESDD, (F77_CONST_CHAR_ARG2 (&jobz, 1), m, n, tmp_data, m1, s_vec, u, m1, vt, nrow_vt1, work.fortran_vec (), lwork, rwork.fortran_vec (), iwork, info F77_CHAR_ARG_LEN (1))); } else assert (0); // impossible if (! (jobv == 'N' || jobv == 'O')) right_sm = right_sm.hermitian (); return info; }