CPPLapack
 All Classes Files Functions Variables Friends
dgematrix-lapack.hpp
Go to the documentation of this file.
00001 //=============================================================================
00002 /*! solve A*X=Y using dgesv\n
00003   The argument is dgematrix Y. Y is overwritten and become the solution X.
00004   A is also overwritten and become P*l*U.
00005 */
00006 inline long dgematrix::dgesv(dgematrix& mat)
00007 {VERBOSE_REPORT;
00008 #ifdef  CPPL_DEBUG
00009   if(m!=n || m!=mat.m){
00010     ERROR_REPORT;
00011     std::cerr << "These two matrices cannot be solved." << std::endl
00012               << "Your input was (" << m << "x" << n << ") and (" << mat.m << "x" << mat.n << ")." << std::endl;
00013     exit(1);
00014   }
00015 #endif//CPPL_DEBUG 
00016   
00017   long NRHS(mat.n), LDA(n), *IPIV(new long[n]), LDB(mat.m), INFO(1);
00018   dgesv_(n, NRHS, array, LDA, IPIV, mat.array, LDB, INFO);
00019   delete [] IPIV;
00020   
00021   if(INFO!=0){
00022     WARNING_REPORT;
00023     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00024   }
00025   return INFO;
00026 }
00027 
00028 //=============================================================================
00029 /*! solve A*x=y using dgesv\n
00030   The argument is dcovector y. y is overwritten and become the solution x.
00031   A is also overwritten and become P*l*U.
00032 */
00033 inline long dgematrix::dgesv(dcovector& vec)
00034 {VERBOSE_REPORT;
00035 #ifdef  CPPL_DEBUG
00036   if(m!=n || m!=vec.l){
00037     ERROR_REPORT;
00038     std::cerr << "These matrix and vector cannot be solved." << std::endl
00039               << "Your input was (" << m << "x" << n << ") and (" << vec.l << ")." << std::endl;
00040     exit(1);
00041   }
00042 #endif//CPPL_DEBUG 
00043   
00044   long NRHS(1), LDA(n), *IPIV(new long[n]), LDB(vec.l), INFO(1);
00045   dgesv_(n, NRHS, array, LDA, IPIV, vec.array, LDB, INFO);
00046   delete [] IPIV;
00047   
00048   if(INFO!=0){
00049     WARNING_REPORT;
00050     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00051   }
00052   return INFO;
00053 }
00054 
00055 ///////////////////////////////////////////////////////////////////////////////
00056 ///////////////////////////////////////////////////////////////////////////////
00057 ///////////////////////////////////////////////////////////////////////////////
00058 
00059 //=============================================================================
00060 /*! solve overdetermined or underdetermined A*X=Y using dgels\n
00061  */
00062 inline long dgematrix::dgels(dgematrix& mat)
00063 {VERBOSE_REPORT;
00064 #ifdef  CPPL_DEBUG
00065   if(m!=mat.m){
00066     ERROR_REPORT;
00067     std::cerr << "These two matrices cannot be solved." << std::endl
00068               << "Your input was (" << m << "x" << n << ") and (" << mat.m << "x" << mat.n << ")." << std::endl;
00069     exit(1);
00070   }
00071 #endif//CPPL_DEBUG    
00072   
00073   if(m<n){ //underdetermined
00074     dgematrix tmp(n,mat.n);
00075     for(long i=0; i<mat.m; i++){ for(long j=0; j<mat.n; j++){
00076       tmp(i,j) =mat(i,j);
00077     }}
00078     mat.clear();
00079     swap(mat,tmp);
00080   }
00081   
00082   char TRANS('n');
00083   long NRHS(mat.n), LDA(m), LDB(mat.m),
00084     LWORK(std::min(m,n)+std::max(std::min(m,n),NRHS)), INFO(1);
00085   double *WORK(new double[LWORK]);
00086   dgels_(TRANS, m, n, NRHS, array, LDA, mat.array, LDB, WORK, LWORK, INFO);
00087   delete [] WORK;
00088   
00089   if(m>n){ //overdetermined
00090     dgematrix tmp(n,mat.n);
00091     for(long i=0; i<tmp.m; i++){ for(long j=0; j<tmp.n; j++){
00092       tmp(i,j) =mat(i,j);
00093     }}
00094     mat.clear();
00095     swap(mat,tmp);
00096   }
00097   
00098   if(INFO!=0){
00099     WARNING_REPORT;
00100     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00101   }
00102   return INFO;
00103 }
00104 
00105 //=============================================================================
00106 /*! solve overdetermined or underdetermined A*x=y using dgels\n
00107  */
00108 inline long dgematrix::dgels(dcovector& vec)
00109 {VERBOSE_REPORT;
00110 #ifdef  CPPL_DEBUG
00111   if(m!=vec.l){
00112     ERROR_REPORT;
00113     std::cerr << "These matrix and vector cannot be solved." << std::endl
00114               << "Your input was (" << m << "x" << n << ") and (" << vec.l << ")." << std::endl;
00115     exit(1);
00116   }
00117 #endif//CPPL_DEBUG    
00118   
00119   if(m<n){ //underdetermined
00120     dcovector tmp(n);
00121     for(long i=0; i<vec.l; i++){ tmp(i)=vec(i); }
00122     vec.clear();
00123     swap(vec,tmp);
00124   }
00125   
00126   char TRANS('n');
00127   long NRHS(1), LDA(m), LDB(vec.l),
00128     LWORK(std::min(m,n)+std::max(std::min(m,n),NRHS)), INFO(1);
00129   double *WORK(new double[LWORK]);
00130   dgels_(TRANS, m, n, NRHS, array, LDA, vec.array, LDB, WORK, LWORK, INFO);
00131   delete [] WORK;
00132   
00133   if(m>n){ //overdetermined
00134     dcovector tmp(n);
00135     for(long i=0; i<tmp.l; i++){ tmp(i)=vec(i); }
00136     vec.clear();
00137     swap(vec,tmp);
00138   }
00139   
00140   if(INFO!=0){
00141     WARNING_REPORT;
00142     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00143   }
00144   return INFO;
00145 }
00146 
00147 //=============================================================================
00148 /*! solve overdetermined or underdetermined A*X=Y using dgels
00149   with the sum of residual squares output\n
00150   The residual is set as the columnwise sum of residual squares 
00151   for overdetermined problems
00152   while it is always zero for underdetermined problems.
00153 */
00154 inline long dgematrix::dgels(dgematrix& mat, drovector& residual)
00155 {VERBOSE_REPORT;
00156 #ifdef  CPPL_DEBUG
00157   if(m!=mat.m){
00158     ERROR_REPORT;
00159     std::cerr << "These two matrices cannot be solved." << std::endl
00160               << "Your input was (" << m << "x" << n << ") and (" << mat.m << "x" << mat.n << ")." << std::endl;
00161     exit(1);
00162   }
00163 #endif//CPPL_DEBUG
00164   
00165   residual.resize(mat.n); residual.zero();
00166   
00167   if(m<n){ //underdetermined
00168     dgematrix tmp(n,mat.n);
00169     for(long i=0; i<mat.m; i++){ for(long j=0; j<mat.n; j++){
00170       tmp(i,j) =mat(i,j);
00171     }}
00172     mat.clear();
00173     swap(mat,tmp);
00174   }
00175   
00176   char TRANS('n');
00177   long NRHS(mat.n), LDA(m), LDB(mat.m),
00178     LWORK(std::min(m,n)+std::max(std::min(m,n),NRHS)), INFO(1);
00179   double *WORK(new double[LWORK]);
00180   dgels_(TRANS, m, n, NRHS, array, LDA, mat.array, LDB, WORK, LWORK, INFO);
00181   delete [] WORK;
00182   
00183   if(m>n){ //overdetermined
00184     for(long i=0; i<residual.l; i++){ for(long j=0; j<m-n; j++){
00185       residual(i) += std::pow(mat(n+j,i), 2.0);
00186     }}
00187     
00188     dgematrix tmp(n,mat.n);
00189     for(long i=0; i<tmp.m; i++){ for(long j=0; j<tmp.n; j++){
00190       tmp(i,j) =mat(i,j);
00191     }}
00192     mat.clear();
00193     swap(mat,tmp);
00194   }
00195   
00196   if(INFO!=0){
00197     WARNING_REPORT;
00198     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00199   }
00200   return INFO;
00201 }
00202 
00203 //=============================================================================
00204 /*! solve overdetermined or underdetermined A*x=y using dgels
00205   with the sum of residual squares output\n
00206   The residual is set as the sum of residual squares 
00207   for overdetermined problems
00208   while it is always zero for underdetermined problems.
00209 */
00210 inline long dgematrix::dgels(dcovector& vec, double& residual)
00211 {VERBOSE_REPORT;
00212 #ifdef  CPPL_DEBUG
00213   if(m!=vec.l){
00214     ERROR_REPORT;
00215     std::cerr << "These matrix and vector cannot be solved." << std::endl
00216               << "Your input was (" << m << "x" << n << ") and (" << vec.l << ")." << std::endl;
00217     exit(1);
00218   }
00219 #endif//CPPL_DEBUG    
00220   
00221   residual=0.0;
00222   
00223   if(m<n){ //underdetermined
00224     dcovector tmp(n);
00225     for(long i=0; i<vec.l; i++){ tmp(i)=vec(i); }
00226     vec.clear();
00227     swap(vec,tmp);
00228   }
00229   
00230   char TRANS('n');
00231   long NRHS(1), LDA(m), LDB(vec.l),
00232     LWORK(std::min(m,n)+std::max(std::min(m,n),NRHS)), INFO(1);
00233   double *WORK(new double[LWORK]);
00234   dgels_(TRANS, m, n, NRHS, array, LDA, vec.array, LDB, WORK, LWORK, INFO);
00235   delete [] WORK;
00236   
00237   if(m>n){ //overdetermined
00238     for(long i=0; i<m-n; i++){ residual+=std::pow(vec(n+i),2.0); }
00239     
00240     dcovector tmp(n);
00241     for(long i=0; i<tmp.l; i++){ tmp(i)=vec(i); }
00242     vec.clear();
00243     swap(vec,tmp);
00244   }
00245   
00246   if(INFO!=0){
00247     WARNING_REPORT;
00248     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00249   }
00250   return INFO;
00251 }
00252 
00253 ///////////////////////////////////////////////////////////////////////////////
00254 ///////////////////////////////////////////////////////////////////////////////
00255 ///////////////////////////////////////////////////////////////////////////////
00256 
00257 //=============================================================================
00258 /*! calculate the least-squares-least-norm solution for overdetermined or 
00259   underdetermined A*x=y using dgelss\n
00260 */
00261 inline long dgematrix::dgelss(dcovector& B, dcovector& S, long& RANK,
00262                               const double RCOND =-1. )
00263 {VERBOSE_REPORT;
00264 #ifdef  CPPL_DEBUG
00265   if(m!=B.l){
00266     ERROR_REPORT;
00267     std::cerr << "These matrix and vector cannot be solved." << std::endl
00268               << "Your input was (" << m << "x" << n << ") and (" << B.l << ")." << std::endl;
00269     exit(1);
00270   }
00271 #endif//CPPL_DEBUG    
00272   
00273   if(m<n){ //underdetermined
00274     dcovector tmp(n);
00275     for(long i=0; i<B.l; i++){ tmp(i)=B(i); }
00276     B.clear();
00277     swap(B,tmp);
00278   }
00279   
00280   S.resize(std::min(m,n));
00281   
00282   long NRHS(1), LDA(m), LDB(B.l),
00283     LWORK(3*std::min(m,n)+std::max(std::max(2*std::min(m,n),std::max(m,n)), NRHS)), INFO(1);
00284   double *WORK(new double[LWORK]);
00285   dgelss_(m, n, NRHS, array, LDA, B.array, LDB, S.array,
00286           RCOND, RANK, WORK, LWORK, INFO);
00287   delete [] WORK;
00288 
00289   if(INFO!=0){
00290     WARNING_REPORT;
00291     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00292   }
00293   return INFO;
00294 }
00295 
00296 //=============================================================================
00297 /*! calculate the least-squares-least-norm solution for overdetermined or 
00298   underdetermined A*x=y using dgelss\n
00299 */
00300 inline long dgematrix::dgelss(dgematrix& B, dcovector& S, long& RANK,
00301                               const double RCOND =-1. )
00302 {VERBOSE_REPORT;
00303 #ifdef  CPPL_DEBUG
00304   if(m!=B.m){
00305     ERROR_REPORT;
00306     std::cerr << "These matrix and vector cannot be solved." << std::endl
00307               << "Your input was (" << m << "x" << n << ") and (" << B.m << "x" << B.n  << ")." << std::endl;
00308     exit(1);
00309   }
00310 #endif//CPPL_DEBUG    
00311   
00312   if(m<n){ //underdetermined
00313     dgematrix tmp(n,B.n);
00314     for(long i=0; i<B.m; i++){
00315       for(long j=0; j<B.n; j++){
00316         tmp(i,j)=B(i,j);
00317       }
00318     }
00319     B.clear();
00320     swap(B,tmp);
00321   }
00322   
00323   S.resize(std::min(m,n));
00324   
00325   long NRHS(B.n), LDA(m), LDB(B.m),
00326     LWORK(3*std::min(m,n)+std::max(std::max(2*std::min(m,n),std::max(m,n)), NRHS)), INFO(1);
00327   double *WORK(new double[LWORK]);
00328   dgelss_(m, n, NRHS, array, LDA, B.array, LDB, S.array,
00329           RCOND, RANK, WORK, LWORK, INFO);
00330   delete [] WORK;
00331 
00332   if(INFO!=0){
00333     WARNING_REPORT;
00334     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00335   }
00336   return INFO;
00337 }
00338 
00339 ///////////////////////////////////////////////////////////////////////////////
00340 ///////////////////////////////////////////////////////////////////////////////
00341 ///////////////////////////////////////////////////////////////////////////////
00342 
00343 //=============================================================================
00344 /*! calculate eigenvalues\n
00345   All of the arguments need not to be initialized. 
00346   wr and wi are overwitten and become 
00347   real and imaginary part of eigenvalues, respectively. 
00348   This matrix is also overwritten. 
00349 */
00350 inline long dgematrix::dgeev(std::vector<double>& wr, std::vector<double>& wi)
00351 {VERBOSE_REPORT;
00352 #ifdef  CPPL_DEBUG
00353   if(m!=n){
00354     ERROR_REPORT;
00355     std::cerr << "This matrix is not a square matrix." << std::endl
00356               << "This matrix is (" << m << "x" << n << ")." << std::endl;
00357     exit(1);
00358   }
00359 #endif//CPPL_DEBUG
00360   
00361   wr.resize(n); wi.resize(n);
00362   char JOBVL('n'), JOBVR('n');
00363   long LDA(n), LDVL(1), LDVR(1), LWORK(3*n), INFO(1);
00364   double *VL(NULL), *VR(NULL), *WORK(new double[LWORK]);
00365   dgeev_(JOBVL, JOBVR, n, array, LDA, &wr[0], &wi[0], 
00366          VL, LDVL, VR, LDVR, WORK, LWORK, INFO);
00367   delete [] WORK; delete [] VL; delete [] VL;
00368   
00369   if(INFO!=0){
00370     WARNING_REPORT;
00371     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00372   }
00373   return INFO;
00374 }
00375 
00376 //=============================================================================
00377 /*! calculate eigenvalues\n
00378   All of the arguments need not to be initialized. 
00379   w are overwitten and become eigenvalues, respectively.
00380   This matrix is also overwritten. 
00381 */
00382 inline long dgematrix::dgeev(zcovector& w)
00383 {VERBOSE_REPORT;
00384   //// call dgeev ////
00385   std::vector<double> wr, wi;
00386   long INFO =dgeev(wr,wi);
00387   
00388   //// assign ////
00389   w.resize(n);
00390   for(long i=0; i<n; i++){
00391     w(i) =comple(wr[i],wi[i]);
00392   }
00393   
00394   return INFO;
00395 }
00396 
00397 //=============================================================================
00398 /*! calculate right eigenvalues and right eigenvectors\n
00399   All of the arguments need not to be initialized. 
00400   wr, wi, vrr, vri are overwitten and become 
00401   real and imaginary part of right eigenvalues and right eigenvectors, 
00402   respectively. 
00403   This matrix is also overwritten. 
00404 */
00405 inline long dgematrix::dgeev
00406 (
00407  std::vector<double>& wr,
00408  std::vector<double>& wi, 
00409  std::vector<dcovector>& vrr,
00410  std::vector<dcovector>& vri
00411  )
00412 {VERBOSE_REPORT;
00413 #ifdef  CPPL_DEBUG
00414   if(m!=n){
00415     ERROR_REPORT;
00416     std::cerr << "This matrix is not a square matrix." << std::endl
00417               << "This matrix is (" << m << "x" << n << ")." << std::endl;
00418     exit(1);
00419   }
00420 #endif//CPPL_DEBUG
00421   
00422   wr.resize(n); wi.resize(n); vrr.resize(n); vri.resize(n);
00423   for(long i=0; i<n; i++){ vrr[i].resize(n); vri[i].resize(n); }
00424   dgematrix VR(n,n);
00425   char JOBVL('n'), JOBVR('V');
00426   long LDA(n), LDVL(1), LDVR(n), LWORK(4*n), INFO(1);
00427   double *VL(NULL), *WORK(new double[LWORK]);
00428   dgeev_(JOBVL, JOBVR, n, array, LDA, &wr[0], &wi[0], 
00429          VL, LDVL, VR.array, LDVR, WORK, LWORK, INFO);
00430   delete [] WORK; delete [] VL;
00431   
00432   //// forming ////
00433   for(long j=0; j<n; j++){
00434     if(fabs(wi[j])<DBL_MIN){
00435       for(long i=0; i<n; i++){
00436         vrr[j](i) = VR(i,j);
00437         vri[j](i) = 0.0;
00438       }
00439     }
00440     else{
00441       for(long i=0; i<n; i++){
00442         vrr[j](i)   = VR(i,j);
00443         vri[j](i)   = VR(i,j+1);
00444         vrr[j+1](i) = VR(i,j);
00445         vri[j+1](i) =-VR(i,j+1);
00446       }
00447       j++;
00448     }
00449   }
00450   
00451   if(INFO!=0){
00452     WARNING_REPORT;
00453     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00454   }
00455   return INFO;
00456 }
00457 
00458 //=============================================================================
00459 /*! calculate left eigenvalues and left eigenvectors\n
00460   All of the arguments need not to be initialized. 
00461   wr, wi, vrr, vri are overwitten and become 
00462   real and imaginary part of left eigenvalues and left eigenvectors, 
00463   respectively. 
00464   This matrix is also overwritten. 
00465 */
00466 inline long dgematrix::dgeev(std::vector<double>& wr, std::vector<double>& wi, 
00467                              std::vector<drovector>& vlr, 
00468                              std::vector<drovector>& vli)
00469 {VERBOSE_REPORT;
00470 #ifdef  CPPL_DEBUG
00471   if(m!=n){
00472     ERROR_REPORT;
00473     std::cerr << "This matrix is not a square matrix." << std::endl
00474               << "This matrix is (" << m << "x" << n << ")." << std::endl;
00475     exit(1);
00476   }
00477 #endif//CPPL_DEBUG
00478   
00479   wr.resize(n); wi.resize(n); vlr.resize(n); vli.resize(n);
00480   for(long i=0; i<n; i++){ vlr[i].resize(n); vli[i].resize(n); }
00481   dgematrix VL(n,n);
00482   char JOBVL('V'), JOBVR('n');
00483   long LDA(n), LDVL(n), LDVR(1), LWORK(4*n), INFO(1);
00484   double *VR(NULL), *WORK(new double[LWORK]);
00485   dgeev_(JOBVL, JOBVR, n, array, LDA, &wr[0], &wi[0], 
00486          VL.array, LDVL, VR, LDVR, WORK, LWORK, INFO);
00487   delete [] WORK; delete [] VR;
00488 
00489   //// forming ////
00490   for(long j=0; j<n; j++){
00491     if(fabs(wi[j])<DBL_MIN){
00492       for(long i=0; i<n; i++){
00493         vlr[j](i) = VL(i,j);
00494         vli[j](i) = 0.0;
00495       }
00496     }
00497     else{
00498       for(long i=0; i<n; i++){
00499         vlr[j](i)   = VL(i,j);
00500         vli[j](i)   =-VL(i,j+1);
00501         vlr[j+1](i) = VL(i,j);
00502         vli[j+1](i) = VL(i,j+1);
00503       }
00504       j++;
00505     }
00506   }
00507   
00508 
00509   if(INFO!=0){
00510     WARNING_REPORT;
00511     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00512   }
00513   return INFO;
00514 }
00515 
00516 
00517 
00518 ///////////////////////////////////////////////////////////////////////////////
00519 ///////////////////////////////////////////////////////////////////////////////
00520 ///////////////////////////////////////////////////////////////////////////////
00521 
00522 //=============================================================================
00523 /*! calculate generalized eigenvalues\n
00524   All of the arguments don't need to be initialized. 
00525   wr and wi are overwitten and become 
00526   real and imaginary part of generalized eigenvalues, respectively. 
00527   This matrix and matB are also overwritten. 
00528 */
00529 inline long dgematrix::dggev(dgematrix& matB,
00530                              std::vector<double>& wr, std::vector<double>& wi)
00531 {VERBOSE_REPORT;
00532 #ifdef  CPPL_DEBUG
00533   if(m!=n){
00534     ERROR_REPORT;
00535     std::cerr << "This matrix is not a square matrix." << std::endl
00536               << "This matrix is (" << m << "x" << n << ")." << std::endl;
00537     exit(1);
00538   }
00539   if(matB.m!=n || matB.n!=n){
00540     ERROR_REPORT;
00541     std::cerr << "The matrix B is not a square matrix having the same size as \"this\" matrix." << std::endl
00542               << "The B matrix is (" << matB.m << "x" << matB.n << ")." << std::endl;
00543     exit(1);
00544   }
00545 #endif//CPPL_DEBUG
00546   
00547   wr.resize(n); wi.resize(n);
00548   char JOBVL('n'), JOBVR('n');
00549   long LDA(n), LDB(n), LDVL(1), LDVR(1), LWORK(8*n), INFO(1);
00550   double *BETA(new double[n]), *VL(NULL), *VR(NULL), 
00551     *WORK(new double[LWORK]);
00552   dggev_(JOBVL, JOBVR, n, array, LDA, matB.array, LDB, &wr[0], &wi[0], 
00553          BETA, VL, LDVL, VR, LDVR, WORK, LWORK, INFO);
00554   delete [] WORK; delete [] VL; delete [] VR;
00555   
00556   //// reforming ////
00557   for(long i=0; i<n; i++){ wr[i]/=BETA[i];  wi[i]/=BETA[i]; }
00558   delete [] BETA; 
00559   
00560   if(INFO!=0){
00561     WARNING_REPORT;
00562     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00563   }
00564   return INFO;
00565 }
00566 
00567 //=============================================================================
00568 /*! calculate generalized eigenvalues and generalized right eigenvectors\n
00569   All of the arguments don't need to be initialized.
00570   wr, wi, vrr and vri are overwitten and become 
00571   real and imaginary part of generalized eigenvalue 
00572   and generalized right eigenvector, respectively. 
00573   This matrix and matB are also overwritten.
00574 */
00575 inline long dgematrix::dggev(dgematrix& matB, 
00576                              std::vector<double>& wr, std::vector<double>& wi, 
00577                              std::vector<dcovector>& vrr, 
00578                              std::vector<dcovector>& vri)
00579 {VERBOSE_REPORT;
00580 #ifdef  CPPL_DEBUG
00581   if(m!=n){
00582     ERROR_REPORT;
00583     std::cerr << "This matrix is not a square matrix." << std::endl
00584               << "This matrix is (" << m << "x" << n << ")." << std::endl;
00585     exit(1);
00586   }
00587   if(matB.m!=n || matB.n!=n){
00588     ERROR_REPORT;
00589     std::cerr << "The matrix B is not a square matrix having the same size as \"this\" matrix." << std::endl
00590               << "The B matrix is (" << matB.m << "x" << matB.n << ")." << std::endl;
00591     exit(1);
00592   }
00593 #endif//CPPL_DEBUG
00594   
00595   wr.resize(n); wi.resize(n); vrr.resize(n); vri.resize(n);
00596   for(long i=0; i<n; i++){ vrr[i].resize(n); vri[i].resize(n); }
00597   dgematrix VR(n,n);
00598   char JOBVL('n'), JOBVR('V');
00599   long LDA(n), LDB(n), LDVL(1), LDVR(n), LWORK(8*n), INFO(1);
00600   double *BETA(new double[n]), *VL(NULL), *WORK(new double[LWORK]);
00601   dggev_(JOBVL, JOBVR, n, array, LDA, matB.array, LDB, &wr[0], &wi[0], 
00602          BETA, VL, LDVL, VR.array, LDVR, WORK, LWORK, INFO);
00603   delete [] WORK; delete [] VL;
00604   
00605   //// reforming ////
00606   for(long i=0; i<n; i++){ wr[i]/=BETA[i];  wi[i]/=BETA[i]; }
00607   delete [] BETA; 
00608   
00609   //// forming ////
00610   for(long j=0; j<n; j++){
00611     if(fabs(wi[j])<DBL_MIN){
00612       for(long i=0; i<n; i++){
00613         vrr[j](i) = VR(i,j);
00614         vri[j](i) = 0.0;
00615       }
00616     }
00617     else{
00618       for(long i=0; i<n; i++){
00619         vrr[j](i)   = VR(i,j);
00620         vri[j](i)   = VR(i,j+1);
00621         vrr[j+1](i) = VR(i,j);
00622         vri[j+1](i) =-VR(i,j+1);
00623       }
00624       j++;
00625     }
00626   }
00627   
00628   if(INFO!=0){
00629     WARNING_REPORT;
00630     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00631   }
00632   return INFO;
00633 }
00634 
00635 //=============================================================================
00636 /*! calculate generalized eigenvalues and generalized left eigenvectors\n
00637   All of the arguments don't need to be initialized.
00638   wr, wi, vlr and vli are overwitten and become 
00639   real and imaginary part of generalized eigenvalue 
00640   and generalized left eigenvector, respectively. 
00641   This matrix and matB are also overwritten.
00642 */
00643 inline long dgematrix::dggev(dgematrix& matB, 
00644                              std::vector<double>& wr, std::vector<double>& wi, 
00645                              std::vector<drovector>& vlr, 
00646                              std::vector<drovector>& vli)
00647 {VERBOSE_REPORT;
00648 #ifdef  CPPL_DEBUG
00649   if(m!=n){
00650     ERROR_REPORT;
00651     std::cerr << "This matrix is not a square matrix." << std::endl
00652               << "This matrix is (" << m << "x" << n << ")." << std::endl;
00653     exit(1);
00654   }
00655   if(matB.m!=n || matB.n!=n){
00656     ERROR_REPORT;
00657     std::cerr << "The matrix B is not a square matrix having the same size as \"this\" matrix." << std::endl
00658               << "The B matrix is (" << matB.m << "x" << matB.n << ")." << std::endl;
00659     exit(1);
00660   }
00661 #endif//CPPL_DEBUG
00662   
00663   wr.resize(n); wi.resize(n); vlr.resize(n); vli.resize(n);
00664   for(long i=0; i<n; i++){ vlr[i].resize(n); vli[i].resize(n); }
00665   dgematrix VL(n,n);
00666   char JOBVL('V'), JOBVR('n');
00667   long LDA(n), LDB(n), LDVL(n), LDVR(1), LWORK(8*n), INFO(1);
00668   double *BETA(new double[n]), *VR(NULL), *WORK(new double[LWORK]);
00669   dggev_(JOBVL, JOBVR, n, array, LDA, matB.array, LDB, &wr[0], &wi[0], 
00670          BETA, VL.array, LDVL, VR, LDVR, WORK, LWORK, INFO);
00671   delete [] WORK; delete [] VR;
00672   
00673   //// reforming ////
00674   for(long i=0; i<n; i++){ wr[i]/=BETA[i];  wi[i]/=BETA[i]; }
00675   delete [] BETA; 
00676   
00677   //// forming ////
00678   for(long j=0; j<n; j++){
00679     if(fabs(wi[j])<DBL_MIN){
00680       for(long i=0; i<n; i++){
00681         vlr[j](i) = VL(i,j);
00682         vli[j](i) = 0.0;
00683       }
00684     }
00685     else{
00686       for(long i=0; i<n; i++){
00687         vlr[j](i)   = VL(i,j);
00688         vli[j](i)   =-VL(i,j+1);
00689         vlr[j+1](i) = VL(i,j);
00690         vli[j+1](i) = VL(i,j+1);
00691       }
00692       j++;
00693     }
00694   }
00695   
00696   if(INFO!=0){
00697     WARNING_REPORT;
00698     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00699   }
00700   return INFO;
00701 }
00702 
00703 ///////////////////////////////////////////////////////////////////////////////
00704 ///////////////////////////////////////////////////////////////////////////////
00705 ///////////////////////////////////////////////////////////////////////////////
00706 
00707 //=============================================================================
00708 /*! compute the singular value decomposition (SVD)\n
00709   The argument is dgbmatrix S.
00710   S doesn't need to be initialized.
00711   S is overwitten and become singular values.
00712   This matrix also overwritten.
00713 */
00714 inline long dgematrix::dgesvd
00715 (
00716  dgbmatrix& S
00717  )
00718 {VERBOSE_REPORT;
00719   char JOBU('n'), JOBVT('n');
00720   long LDA(m), LDU(m), LDVT(n), LWORK(std::max(3*std::min(m,n)+std::max(m,n),5*std::min(m,n))), INFO(1);
00721   double *WORK(new double[LWORK]);
00722   S.resize(m,n,0,0);
00723   
00724   dgesvd_(JOBU, JOBVT, m, n, array, LDA, S.array, NULL, LDU, NULL, LDVT, WORK, LWORK, INFO);
00725   delete [] WORK;
00726   
00727   if(INFO!=0){
00728     WARNING_REPORT;
00729     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00730   }
00731   return INFO;
00732 }
00733 
00734 //=============================================================================
00735 /*! compute the singular value decomposition (SVD)\n
00736   The arguments are dcocector S, dgematrix U and VT.
00737   All of them need not to be initialized.
00738   S, U and VT are overwitten and become singular values, 
00739   left singular vectors,
00740   and right singular vectors respectively.
00741   This matrix also overwritten.
00742 */
00743 inline long dgematrix::dgesvd(dcovector& S, dgematrix& U, dgematrix& VT)
00744 {VERBOSE_REPORT;
00745   char JOBU('A'), JOBVT('A');
00746   long LDA(m), LDU(m), LDVT(n),
00747     LWORK(std::max(3*std::min(m,n)+std::max(m,n),5*std::min(m,n))), INFO(1);
00748   double *WORK(new double[LWORK]);
00749   S.resize(std::min(m,n)); U.resize(LDU,m); VT.resize(LDVT,n);
00750   
00751   dgesvd_(JOBU, JOBVT, m, n, array, LDA, S.array, U.array,
00752           LDU, VT.array, LDVT, WORK, LWORK, INFO);
00753   delete [] WORK;
00754   
00755   if(INFO!=0){
00756     WARNING_REPORT;
00757     std::cerr << "Serious trouble happend. INFO = " << INFO << "." << std::endl;
00758   }
00759   return INFO;
00760 }
00761 
00762 ///////////////////////////////////////////////////////////////////////////////
00763 ///////////////////////////////////////////////////////////////////////////////
00764 ///////////////////////////////////////////////////////////////////////////////
00765 
00766 //=============================================================================
00767 /*! solve the linear equality-constrained least squares (LSE) problem\n
00768   Input matrix and vectors, B, c, and d, are overwitten.
00769   This matrix is also overwritten.
00770   The solution vector x is to be automatically resized.
00771 */
00772 inline long dgematrix::dgglse
00773 (
00774  dgematrix& B,
00775  dcovector& c,
00776  dcovector& d,
00777  dcovector& x
00778  )
00779 {VERBOSE_REPORT;
00780 #ifdef  CPPL_DEBUG
00781   if(m!=c.l){
00782     ERROR_REPORT;
00783     std::cerr << "A.m and c.l should be the same." << std::endl
00784               << "Your input was A.m=" << m << " and c.l=" << c.l << std::endl;
00785     exit(1);
00786   }
00787   if(B.m!=d.l){
00788     ERROR_REPORT;
00789     std::cerr << "B.m and d.l should be the same." << std::endl
00790               << "Your input was B.m=" << B.m << " and d.l=" << d.l << std::endl;
00791     exit(1);
00792   }
00793   if( !(B.m<=n) || !(n<=m+B.m) ){
00794     ERROR_REPORT;
00795     std::cerr << "B.m<=A.n<=A.m+B.m should be satisfied." << std::endl
00796               << "Your input was B.m=" << B.m << ", A.n=" << n << ", and A.m+B.m=" << m+B.m << std::endl;
00797     exit(1);
00798   }
00799 #endif//CPPL_DEBUG
00800   
00801   long lwork(-1), info(1);
00802   dcovector work(1);
00803   x.resize(n);
00804   
00805   //////// workspace query ////////
00806   dgglse_(m, n, B.m, array, std::max(1l,m), B.array, std::max(1l,B.m), c.array, d.array, x.array, work.array, lwork, info);
00807   lwork =long(work(0));
00808   work.resize(lwork);
00809   info =1;
00810   
00811   //////// solve ////////
00812   dgglse_(m, n, B.m, array, std::max(1l,m), B.array, std::max(1l,B.m), c.array, d.array, x.array, work.array, lwork,  info);
00813   work.clear();
00814   
00815   if(info!=0){
00816     WARNING_REPORT;
00817     std::cerr << "Serious trouble happend. info = " << info << "." << std::endl;
00818   }
00819   return info;
00820 }
 All Classes Files Functions Variables Friends