| [0b990d] | 1 | //
 | 
|---|
 | 2 | // svd.cc
 | 
|---|
 | 3 | //
 | 
|---|
 | 4 | // Copyright (C) 2005 Edward Valeev
 | 
|---|
 | 5 | //
 | 
|---|
 | 6 | // Author: Edward Valeev <edward.valeev@chemistry.gatech.edu>
 | 
|---|
 | 7 | // Maintainer: EV
 | 
|---|
 | 8 | //
 | 
|---|
 | 9 | // This file is part of the SC Toolkit.
 | 
|---|
 | 10 | //
 | 
|---|
 | 11 | // The SC Toolkit is free software; you can redistribute it and/or modify
 | 
|---|
 | 12 | // it under the terms of the GNU Library General Public License as published by
 | 
|---|
 | 13 | // the Free Software Foundation; either version 2, or (at your option)
 | 
|---|
 | 14 | // any later version.
 | 
|---|
 | 15 | //
 | 
|---|
 | 16 | // The SC Toolkit is distributed in the hope that it will be useful,
 | 
|---|
 | 17 | // but WITHOUT ANY WARRANTY; without even the implied warranty of
 | 
|---|
 | 18 | // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 | 
|---|
 | 19 | // GNU Library General Public License for more details.
 | 
|---|
 | 20 | //
 | 
|---|
 | 21 | // You should have received a copy of the GNU Library General Public License
 | 
|---|
 | 22 | // along with the SC Toolkit; see the file COPYING.LIB.  If not, write to
 | 
|---|
 | 23 | // the Free Software Foundation, 675 Mass Ave, Cambridge, MA 02139, USA.
 | 
|---|
 | 24 | //
 | 
|---|
 | 25 | // The U.S. Government is granted a limited license as per AL 91-7.
 | 
|---|
 | 26 | //
 | 
|---|
 | 27 | 
 | 
|---|
 | 28 | #ifdef __GNUG__
 | 
|---|
 | 29 | #pragma implementation
 | 
|---|
 | 30 | #endif
 | 
|---|
 | 31 | 
 | 
|---|
 | 32 | #include <cmath>
 | 
|---|
 | 33 | #include <stdexcept>
 | 
|---|
 | 34 | #include <util/misc/formio.h>
 | 
|---|
 | 35 | #include <chemistry/qc/mbptr12/lapack.h>
 | 
|---|
 | 36 | #include <chemistry/qc/mbptr12/svd.h>
 | 
|---|
 | 37 | 
 | 
|---|
 | 38 | using namespace std;
 | 
|---|
 | 39 | using namespace sc;
 | 
|---|
 | 40 | 
 | 
|---|
 | 41 | namespace sc {
 | 
|---|
 | 42 | 
 | 
|---|
 | 43 |   namespace exp {
 | 
|---|
 | 44 | 
 | 
|---|
 | 45 |     /** Uses LAPACK's DGESVD to perform SVD of A:
 | 
|---|
 | 46 |         A = U * Sigma * V
 | 
|---|
 | 47 |     */
 | 
|---|
 | 48 |     void lapack_svd(const RefSCMatrix& A,
 | 
|---|
 | 49 |                     RefSCMatrix& U,
 | 
|---|
 | 50 |                     RefDiagSCMatrix& Sigma,
 | 
|---|
 | 51 |                     RefSCMatrix& V)
 | 
|---|
 | 52 |     {
 | 
|---|
 | 53 |       /*
 | 
|---|
 | 54 |         In terms of C matrixes, Fortran77 call is to compute SVD of A^t
 | 
|---|
 | 55 | 
 | 
|---|
 | 56 |           A^t = V^t * Sigma^t * U^t
 | 
|---|
 | 57 | 
 | 
|---|
 | 58 |         DGESVD uses notation A = U * Sigma * Vt, hence the C-F77 correspondence is as follows:
 | 
|---|
 | 59 | 
 | 
|---|
 | 60 |         C     F77
 | 
|---|
 | 61 |        ---    ---
 | 
|---|
 | 62 |         A      A^t
 | 
|---|
 | 63 |         V      U
 | 
|---|
 | 64 |         Sigma  Sigma
 | 
|---|
 | 65 |         U      Vt
 | 
|---|
 | 66 | 
 | 
|---|
 | 67 |       */
 | 
|---|
 | 68 |       int m = A.nrow();
 | 
|---|
 | 69 |       int n = A.ncol();
 | 
|---|
 | 70 |       int nsigma = m<n ? m : n;
 | 
|---|
 | 71 |       int lwork = m<n ? 5*n : 5*m;
 | 
|---|
 | 72 |       char jobvt = U ? 'A' : 'N';
 | 
|---|
 | 73 |       char jobu = V ? 'A' : 'N';
 | 
|---|
 | 74 | 
 | 
|---|
 | 75 |       double* a_data = new double[m*n];
 | 
|---|
 | 76 |       A.convert(a_data);
 | 
|---|
 | 77 | 
 | 
|---|
 | 78 |       double* u_data = 0;
 | 
|---|
 | 79 |       if (U) {
 | 
|---|
 | 80 |         u_data = new double[m*m];
 | 
|---|
 | 81 |       }
 | 
|---|
 | 82 |       double* v_data = 0;
 | 
|---|
 | 83 |       if (V) {
 | 
|---|
 | 84 |         v_data = new double[n*n];
 | 
|---|
 | 85 |       }
 | 
|---|
 | 86 |       double* s_data = new double[nsigma];
 | 
|---|
 | 87 |       double *work = new double[lwork];
 | 
|---|
 | 88 |       int info = 0;
 | 
|---|
 | 89 | 
 | 
|---|
 | 90 |       F77_DGESVD(&jobu, &jobvt, &n, &m,
 | 
|---|
 | 91 |                  a_data, &n, s_data,
 | 
|---|
 | 92 |                  v_data, &n,
 | 
|---|
 | 93 |                  u_data, &m,
 | 
|---|
 | 94 |                  work, &lwork, &info);
 | 
|---|
 | 95 | 
 | 
|---|
 | 96 |       if (info !=0) {
 | 
|---|
 | 97 |         delete[] work;
 | 
|---|
 | 98 |         delete[] a_data;
 | 
|---|
 | 99 |         delete[] s_data;
 | 
|---|
 | 100 |         if (u_data) delete[] u_data;
 | 
|---|
 | 101 |         if (v_data) delete[] v_data;
 | 
|---|
 | 102 |         throw std::runtime_error("lapack_svd -- call to LAPACK routine failed");
 | 
|---|
 | 103 |       }
 | 
|---|
 | 104 | 
 | 
|---|
 | 105 |       Sigma.assign(s_data);
 | 
|---|
 | 106 |       if (U) {
 | 
|---|
 | 107 |         U.assign(u_data);
 | 
|---|
 | 108 |       }
 | 
|---|
 | 109 |       if (V) {
 | 
|---|
 | 110 |         V.assign(v_data);
 | 
|---|
 | 111 |       }
 | 
|---|
 | 112 |       
 | 
|---|
 | 113 |       delete[] work;
 | 
|---|
 | 114 |       delete[] a_data;
 | 
|---|
 | 115 |       delete[] s_data;
 | 
|---|
 | 116 |       if (u_data) delete[] u_data;
 | 
|---|
 | 117 |       if (v_data) delete[] v_data;
 | 
|---|
 | 118 |     }
 | 
|---|
 | 119 | 
 | 
|---|
 | 120 | 
 | 
|---|
 | 121 |     /** Uses LAPACK's DSPSVX to solve symmetric non-definite linear system AX = B, where B is a RefSCVector
 | 
|---|
 | 122 |     */
 | 
|---|
 | 123 |     void lapack_linsolv_symmnondef(const RefSymmSCMatrix& A,
 | 
|---|
 | 124 |                                    RefSCVector& X,
 | 
|---|
 | 125 |                                    const RefSCVector& B)
 | 
|---|
 | 126 |     {
 | 
|---|
 | 127 |       int n = A.n();
 | 
|---|
 | 128 |       if (n != B.n())
 | 
|---|
 | 129 |         throw std::runtime_error("lapack_linsolv_symmnondef() -- dimensions of A and B do not match");
 | 
|---|
 | 130 |       if (n != X.n())
 | 
|---|
 | 131 |         throw std::runtime_error("lapack_linsolv_symmnondef() -- dimensions of A and X do not match");
 | 
|---|
 | 132 | 
 | 
|---|
 | 133 |       // convert A to packed upper triangular form
 | 
|---|
 | 134 |       int ntri = n*(n+1)/2;
 | 
|---|
 | 135 |       double* AP = new double[ntri];
 | 
|---|
 | 136 |       A.convert(AP);
 | 
|---|
 | 137 |       double* AFP = new double[ntri];
 | 
|---|
 | 138 |       int* ipiv = new int[n];
 | 
|---|
 | 139 |       double* BB = new double[n];
 | 
|---|
 | 140 |       B.convert(BB);
 | 
|---|
 | 141 |       double* XX = new double[n];
 | 
|---|
 | 142 | 
 | 
|---|
 | 143 |       char fact = 'N';
 | 
|---|
 | 144 |       char uplo = 'U';
 | 
|---|
 | 145 |       int nrhs = 1;
 | 
|---|
 | 146 |       double rcond = 0.0;
 | 
|---|
 | 147 |       double* ferr = new double[1];
 | 
|---|
 | 148 |       double* berr = new double[1];
 | 
|---|
 | 149 |       double* work = new double[3*n];
 | 
|---|
 | 150 |       int* iwork = new int[n];
 | 
|---|
 | 151 |       int info = 0;
 | 
|---|
 | 152 |       
 | 
|---|
 | 153 |       F77_DSPSVX(&fact, &uplo, &n, &nrhs, AP, AFP, ipiv, BB, &n, XX, &n, &rcond, ferr, berr, work, iwork, &info);
 | 
|---|
 | 154 | 
 | 
|---|
 | 155 |       if (info) {
 | 
|---|
 | 156 |         
 | 
|---|
 | 157 |         delete[] ferr;
 | 
|---|
 | 158 |         delete[] berr;
 | 
|---|
 | 159 |         delete[] work;
 | 
|---|
 | 160 |         delete[] iwork;
 | 
|---|
 | 161 | 
 | 
|---|
 | 162 |         delete[] AP;
 | 
|---|
 | 163 |         delete[] AFP;
 | 
|---|
 | 164 |         delete[] BB;
 | 
|---|
 | 165 |         delete[] XX;
 | 
|---|
 | 166 | 
 | 
|---|
 | 167 |         if (info > 0 && info <= n)
 | 
|---|
 | 168 |           throw std::runtime_error("lapack_linsolv_symmnondef() -- matrix A has factors which are exactly zero");
 | 
|---|
 | 169 |         if (info == n + 1)
 | 
|---|
 | 170 |           throw std::runtime_error("lapack_linsolv_symmnondef() -- matrix A is singular");
 | 
|---|
 | 171 |         if (info < 0)
 | 
|---|
 | 172 |           throw std::runtime_error("lapack_linsolv_symmnondef() -- one of the arguments to F77_DSPSVX is invalid");
 | 
|---|
 | 173 |       }
 | 
|---|
 | 174 |       else {
 | 
|---|
 | 175 |         X.assign(XX);
 | 
|---|
 | 176 | 
 | 
|---|
 | 177 |         delete[] ferr;
 | 
|---|
 | 178 |         delete[] berr;
 | 
|---|
 | 179 |         delete[] work;
 | 
|---|
 | 180 |         delete[] iwork;
 | 
|---|
 | 181 | 
 | 
|---|
 | 182 |         delete[] AP;
 | 
|---|
 | 183 |         delete[] AFP;
 | 
|---|
 | 184 |         delete[] BB;
 | 
|---|
 | 185 |         delete[] XX;
 | 
|---|
 | 186 |       }
 | 
|---|
 | 187 |     }
 | 
|---|
 | 188 | 
 | 
|---|
 | 189 |     /** Uses LAPACK's DSPSVX to solve symmetric non-definite linear system AX = B, where B is a RefSCMatrix.
 | 
|---|
 | 190 |         Dimensions of X and B must match.
 | 
|---|
 | 191 |     */
 | 
|---|
 | 192 |     void lapack_linsolv_symmnondef(const RefSymmSCMatrix& A,
 | 
|---|
 | 193 |                                    RefSCMatrix& X,
 | 
|---|
 | 194 |                                    const RefSCMatrix& B)
 | 
|---|
 | 195 |     {
 | 
|---|
 | 196 |       int n = A.n();
 | 
|---|
 | 197 |       int nrhs = B.ncol();
 | 
|---|
 | 198 |       if (n != B.nrow())
 | 
|---|
 | 199 |         throw std::runtime_error("lapack_linsolv_symmnondef() -- dimensions of A and B do not match");
 | 
|---|
 | 200 |       if (n != X.nrow())
 | 
|---|
 | 201 |         throw std::runtime_error("lapack_linsolv_symmnondef() -- dimensions of A and X do not match");
 | 
|---|
 | 202 |       if (X.ncol() != nrhs)
 | 
|---|
 | 203 |         throw std::runtime_error("lapack_linsolv_symmnondef() -- dimensions of B and X do not match");
 | 
|---|
 | 204 | 
 | 
|---|
 | 205 |       // convert A to packed upper triangular form
 | 
|---|
 | 206 |       int ntri = n*(n+1)/2;
 | 
|---|
 | 207 |       double* AP = new double[ntri];
 | 
|---|
 | 208 |       A.convert(AP);
 | 
|---|
 | 209 |       double* AFP = new double[ntri];
 | 
|---|
 | 210 |       int* ipiv = new int[n];
 | 
|---|
 | 211 |       double* BB = new double[nrhs*n];
 | 
|---|
 | 212 |       B.t().convert(BB);
 | 
|---|
 | 213 |       double* XX = new double[nrhs*n];
 | 
|---|
 | 214 | 
 | 
|---|
 | 215 |       char fact = 'N';
 | 
|---|
 | 216 |       char uplo = 'U';
 | 
|---|
 | 217 |       double rcond = 0.0;
 | 
|---|
 | 218 |       double* ferr = new double[nrhs];
 | 
|---|
 | 219 |       double* berr = new double[nrhs];
 | 
|---|
 | 220 |       double* work = new double[3*n];
 | 
|---|
 | 221 |       int* iwork = new int[n];
 | 
|---|
 | 222 |       int info = 0;
 | 
|---|
 | 223 |       
 | 
|---|
 | 224 |       F77_DSPSVX(&fact, &uplo, &n, &nrhs, AP, AFP, ipiv, BB, &n, XX, &n, &rcond, ferr, berr, work, iwork, &info);
 | 
|---|
 | 225 | 
 | 
|---|
 | 226 |       if (info) {
 | 
|---|
 | 227 |         
 | 
|---|
 | 228 |         delete[] ferr;
 | 
|---|
 | 229 |         delete[] berr;
 | 
|---|
 | 230 |         delete[] work;
 | 
|---|
 | 231 |         delete[] iwork;
 | 
|---|
 | 232 | 
 | 
|---|
 | 233 |         delete[] AP;
 | 
|---|
 | 234 |         delete[] AFP;
 | 
|---|
 | 235 |         delete[] BB;
 | 
|---|
 | 236 |         delete[] XX;
 | 
|---|
 | 237 | 
 | 
|---|
 | 238 |         if (info > 0 && info <= n)
 | 
|---|
 | 239 |           throw std::runtime_error("lapack_linsolv_symmnondef() -- matrix A has factors which are exactly zero");
 | 
|---|
 | 240 |         if (info == n + 1)
 | 
|---|
 | 241 |           throw std::runtime_error("lapack_linsolv_symmnondef() -- matrix A is singular");
 | 
|---|
 | 242 |         if (info < 0)
 | 
|---|
 | 243 |           throw std::runtime_error("lapack_linsolv_symmnondef() -- one of the arguments to F77_DSPSVX is invalid");
 | 
|---|
 | 244 |       }
 | 
|---|
 | 245 |       else {
 | 
|---|
 | 246 |         RefSCMatrix Xt = X.t();
 | 
|---|
 | 247 |         Xt.assign(XX);
 | 
|---|
 | 248 |         X = Xt.t();
 | 
|---|
 | 249 |         Xt = 0;
 | 
|---|
 | 250 | 
 | 
|---|
 | 251 |         delete[] ferr;
 | 
|---|
 | 252 |         delete[] berr;
 | 
|---|
 | 253 |         delete[] work;
 | 
|---|
 | 254 |         delete[] iwork;
 | 
|---|
 | 255 | 
 | 
|---|
 | 256 |         delete[] AP;
 | 
|---|
 | 257 |         delete[] AFP;
 | 
|---|
 | 258 |         delete[] BB;
 | 
|---|
 | 259 |         delete[] XX;
 | 
|---|
 | 260 |       }
 | 
|---|
 | 261 |     }
 | 
|---|
 | 262 | 
 | 
|---|
 | 263 |   };
 | 
|---|
 | 264 |   
 | 
|---|
 | 265 | };
 | 
|---|
 | 266 | 
 | 
|---|
 | 267 | /////////////////////////////////////////////////////////////////////////////
 | 
|---|
 | 268 | 
 | 
|---|
 | 269 | // Local Variables:
 | 
|---|
 | 270 | // mode: c++
 | 
|---|
 | 271 | // c-file-style: "CLJ-CONDENSED"
 | 
|---|
 | 272 | // End:
 | 
|---|