sparse-dmsolve.cc

Go to the documentation of this file.
00001 /*
00002 
00003 Copyright (C) 2006-2012 David Bateman
00004 
00005 This file is part of Octave.
00006 
00007 Octave is free software; you can redistribute it and/or modify it
00008 under the terms of the GNU General Public License as published by the
00009 Free Software Foundation; either version 3 of the License, or (at your
00010 option) any later version.
00011 
00012 Octave is distributed in the hope that it will be useful, but WITHOUT
00013 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00014 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00015 for more details.
00016 
00017 You should have received a copy of the GNU General Public License
00018 along with Octave; see the file COPYING.  If not, see
00019 <http://www.gnu.org/licenses/>.
00020 
00021 */
00022 
00023 #ifdef HAVE_CONFIG_H
00024 #include <config.h>
00025 #endif
00026 
00027 #include <vector>
00028 
00029 #include "MArray.h"
00030 #include "MSparse.h"
00031 #include "SparseQR.h"
00032 #include "SparseCmplxQR.h"
00033 #include "MatrixType.h"
00034 #include "oct-sort.h"
00035 #include "oct-locbuf.h"
00036 #include "oct-inttypes.h"
00037 
00038 template <class T>
00039 static MSparse<T>
00040 dmsolve_extract (const MSparse<T> &A, const octave_idx_type *Pinv,
00041                 const octave_idx_type *Q, octave_idx_type rst,
00042                 octave_idx_type rend, octave_idx_type cst,
00043                 octave_idx_type cend, octave_idx_type maxnz = -1,
00044                 bool lazy = false)
00045 {
00046   octave_idx_type nr = rend - rst, nc = cend - cst;
00047   maxnz = (maxnz < 0 ? A.nnz () : maxnz);
00048   octave_idx_type nz;
00049 
00050   // Cast to uint64 to handle overflow in this multiplication
00051   if (octave_uint64 (nr)*octave_uint64 (nc) < octave_uint64 (maxnz))
00052     nz = nr*nc;
00053   else
00054     nz = maxnz;
00055 
00056   MSparse<T> B (nr, nc, (nz < maxnz ? nz : maxnz));
00057   // Some sparse functions can support lazy indexing (where elements
00058   // in the row are in no particular order), even though octave in
00059   // general can't. For those functions that can using it is a big
00060   // win here in terms of speed.
00061   if (lazy)
00062     {
00063       nz = 0;
00064       for (octave_idx_type j = cst ; j < cend ; j++)
00065         {
00066           octave_idx_type qq = (Q ? Q [j] : j);
00067           B.xcidx (j - cst) = nz;
00068           for (octave_idx_type p = A.cidx(qq) ; p < A.cidx (qq+1) ; p++)
00069             {
00070               octave_quit ();
00071               octave_idx_type r = (Pinv ? Pinv [A.ridx (p)] : A.ridx (p));
00072               if (r >= rst && r < rend)
00073                 {
00074                   B.xdata (nz) = A.data (p);
00075                   B.xridx (nz++) =  r - rst ;
00076                 }
00077             }
00078         }
00079       B.xcidx (cend - cst) = nz ;
00080     }
00081   else
00082     {
00083       OCTAVE_LOCAL_BUFFER (T, X, rend - rst);
00084       octave_sort<octave_idx_type> sort;
00085       octave_idx_type *ri = B.xridx();
00086       nz = 0;
00087       for (octave_idx_type j = cst ; j < cend ; j++)
00088         {
00089           octave_idx_type qq = (Q ? Q [j] : j);
00090           B.xcidx (j - cst) = nz;
00091           for (octave_idx_type p = A.cidx(qq) ; p < A.cidx (qq+1) ; p++)
00092             {
00093               octave_quit ();
00094               octave_idx_type r = (Pinv ? Pinv [A.ridx (p)] : A.ridx (p));
00095               if (r >= rst && r < rend)
00096                 {
00097                   X [r-rst] = A.data (p);
00098                   B.xridx (nz++) =  r - rst ;
00099                 }
00100             }
00101           sort.sort (ri + B.xcidx (j - cst), nz - B.xcidx (j - cst));
00102           for (octave_idx_type p = B.cidx (j - cst); p < nz; p++)
00103             B.xdata (p) = X [B.xridx (p)];
00104         }
00105       B.xcidx (cend - cst) = nz ;
00106     }
00107 
00108   return B;
00109 }
00110 
00111 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00112 static MSparse<double>
00113 dmsolve_extract (const MSparse<double> &A, const octave_idx_type *Pinv,
00114                 const octave_idx_type *Q, octave_idx_type rst,
00115                 octave_idx_type rend, octave_idx_type cst,
00116                 octave_idx_type cend, octave_idx_type maxnz,
00117                 bool lazy);
00118 
00119 static MSparse<Complex>
00120 dmsolve_extract (const MSparse<Complex> &A, const octave_idx_type *Pinv,
00121                 const octave_idx_type *Q, octave_idx_type rst,
00122                 octave_idx_type rend, octave_idx_type cst,
00123                 octave_idx_type cend, octave_idx_type maxnz,
00124                 bool lazy);
00125 #endif
00126 
00127 template <class T>
00128 static MArray<T>
00129 dmsolve_extract (const MArray<T> &m, const octave_idx_type *,
00130                  const octave_idx_type *, octave_idx_type r1,
00131                  octave_idx_type r2, octave_idx_type c1,
00132                  octave_idx_type c2)
00133 {
00134   r2 -= 1;
00135   c2 -= 1;
00136   if (r1 > r2) { octave_idx_type tmp = r1; r1 = r2; r2 = tmp; }
00137   if (c1 > c2) { octave_idx_type tmp = c1; c1 = c2; c2 = tmp; }
00138 
00139   octave_idx_type new_r = r2 - r1 + 1;
00140   octave_idx_type new_c = c2 - c1 + 1;
00141 
00142   MArray<T> result (dim_vector (new_r, new_c));
00143 
00144   for (octave_idx_type j = 0; j < new_c; j++)
00145     for (octave_idx_type i = 0; i < new_r; i++)
00146       result.xelem (i, j) = m.elem (r1+i, c1+j);
00147 
00148   return result;
00149 }
00150 
00151 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00152 static MArray<double>
00153 dmsolve_extract (const MArray<double> &m, const octave_idx_type *,
00154                  const octave_idx_type *, octave_idx_type r1,
00155                  octave_idx_type r2, octave_idx_type c1,
00156                  octave_idx_type c2)
00157 
00158 static MArray<Complex>
00159 dmsolve_extract (const MArray<Complex> &m, const octave_idx_type *,
00160                  const octave_idx_type *, octave_idx_type r1,
00161                  octave_idx_type r2, octave_idx_type c1,
00162                  octave_idx_type c2)
00163 #endif
00164 
00165 template <class T>
00166 static void
00167 dmsolve_insert (MArray<T> &a, const MArray<T> &b, const octave_idx_type *Q,
00168                octave_idx_type r, octave_idx_type c)
00169 {
00170   T *ax = a.fortran_vec();
00171   const T *bx = b.fortran_vec();
00172   octave_idx_type anr = a.rows();
00173   octave_idx_type nr = b.rows();
00174   octave_idx_type nc = b.cols();
00175   for (octave_idx_type j = 0; j < nc; j++)
00176     {
00177       octave_idx_type aoff = (c + j) * anr;
00178       octave_idx_type boff = j * nr;
00179       for (octave_idx_type i = 0; i < nr; i++)
00180         {
00181           octave_quit ();
00182           ax [Q [r + i] + aoff] = bx [i + boff];
00183         }
00184     }
00185 }
00186 
00187 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00188 static void
00189 dmsolve_insert (MArray<double> &a, const MArray<double> &b,
00190                const octave_idx_type *Q, octave_idx_type r, octave_idx_type c);
00191 
00192 static void
00193 dmsolve_insert (MArray<Complex> &a, const MArray<Complex> &b,
00194                const octave_idx_type *Q, octave_idx_type r, octave_idx_type c);
00195 #endif
00196 
00197 template <class T>
00198 static void
00199 dmsolve_insert (MSparse<T> &a, const MSparse<T> &b, const octave_idx_type *Q,
00200                octave_idx_type r, octave_idx_type c)
00201 {
00202   octave_idx_type b_rows = b.rows ();
00203   octave_idx_type b_cols = b.cols ();
00204   octave_idx_type nr = a.rows ();
00205   octave_idx_type nc = a.cols ();
00206 
00207   OCTAVE_LOCAL_BUFFER (octave_idx_type, Qinv, nr);
00208   for (octave_idx_type i = 0; i < nr; i++)
00209     Qinv [Q [i]] = i;
00210 
00211   // First count the number of elements in the final array
00212   octave_idx_type nel = a.xcidx(c) + b.nnz ();
00213 
00214   if (c + b_cols < nc)
00215     nel += a.xcidx(nc) - a.xcidx(c + b_cols);
00216 
00217   for (octave_idx_type i = c; i < c + b_cols; i++)
00218     for (octave_idx_type j = a.xcidx(i); j < a.xcidx(i+1); j++)
00219       if (Qinv [a.xridx(j)] < r || Qinv [a.xridx(j)] >= r + b_rows)
00220         nel++;
00221 
00222   OCTAVE_LOCAL_BUFFER (T, X, nr);
00223   octave_sort<octave_idx_type> sort;
00224   MSparse<T> tmp (a);
00225   a = MSparse<T> (nr, nc, nel);
00226   octave_idx_type *ri = a.xridx();
00227 
00228   for (octave_idx_type i = 0; i < tmp.cidx(c); i++)
00229     {
00230       a.xdata(i) = tmp.xdata(i);
00231       a.xridx(i) = tmp.xridx(i);
00232     }
00233   for (octave_idx_type i = 0; i < c + 1; i++)
00234     a.xcidx(i) = tmp.xcidx(i);
00235 
00236   octave_idx_type ii = a.xcidx(c);
00237 
00238   for (octave_idx_type i = c; i < c + b_cols; i++)
00239     {
00240       octave_quit ();
00241 
00242       for (octave_idx_type j = tmp.xcidx(i); j < tmp.xcidx(i+1); j++)
00243         if (Qinv [tmp.xridx(j)] < r ||  Qinv [tmp.xridx(j)] >= r + b_rows)
00244           {
00245             X [tmp.xridx(j)] = tmp.xdata(j);
00246             a.xridx(ii++) = tmp.xridx(j);
00247           }
00248 
00249       octave_quit ();
00250 
00251       for (octave_idx_type j = b.cidx(i-c); j < b.cidx(i-c+1); j++)
00252         {
00253           X [Q [r + b.ridx(j)]] = b.data(j);
00254           a.xridx(ii++) = Q [r + b.ridx(j)];
00255         }
00256 
00257       sort.sort (ri + a.xcidx (i), ii - a.xcidx (i));
00258       for (octave_idx_type p = a.xcidx (i); p < ii; p++)
00259         a.xdata (p) = X [a.xridx (p)];
00260       a.xcidx(i+1) = ii;
00261     }
00262 
00263   for (octave_idx_type i = c + b_cols; i < nc; i++)
00264     {
00265       for (octave_idx_type j = tmp.xcidx(i); j < tmp.cidx(i+1); j++)
00266         {
00267           a.xdata(ii) = tmp.xdata(j);
00268           a.xridx(ii++) = tmp.xridx(j);
00269         }
00270       a.xcidx(i+1) = ii;
00271     }
00272 }
00273 
00274 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00275 static void
00276 dmsolve_insert (MSparse<double> &a, const SparseMatrix &b,
00277                const octave_idx_type *Q, octave_idx_type r, octave_idx_type c);
00278 
00279 static void
00280 dmsolve_insert (MSparse<Complex> &a, const MSparse<Complex> &b,
00281                const octave_idx_type *Q, octave_idx_type r, octave_idx_type c);
00282 #endif
00283 
00284 template <class T, class RT>
00285 static void
00286 dmsolve_permute (MArray<RT> &a, const MArray<T>& b, const octave_idx_type *p)
00287 {
00288   octave_idx_type b_nr = b.rows ();
00289   octave_idx_type b_nc = b.cols ();
00290   const T *Bx = b.fortran_vec();
00291   a.resize (dim_vector (b_nr, b_nc));
00292   RT *Btx = a.fortran_vec();
00293   for (octave_idx_type j = 0; j < b_nc; j++)
00294     {
00295       octave_idx_type off = j * b_nr;
00296       for (octave_idx_type i = 0; i < b_nr; i++)
00297         {
00298           octave_quit ();
00299           Btx [p [i] + off] = Bx [ i + off];
00300         }
00301     }
00302 }
00303 
00304 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00305 static void
00306 dmsolve_permute (MArray<double> &a, const MArray<double>& b,
00307                  const octave_idx_type *p);
00308 
00309 static void
00310 dmsolve_permute (MArray<Complex> &a, const MArray<double>& b,
00311                  const octave_idx_type *p);
00312 
00313 static void
00314 dmsolve_permute (MArray<Complex> &a, const MArray<Complex>& b,
00315                  const octave_idx_type *p);
00316 #endif
00317 
00318 template <class T, class RT>
00319 static void
00320 dmsolve_permute (MSparse<RT> &a, const MSparse<T>& b, const octave_idx_type *p)
00321 {
00322   octave_idx_type b_nr = b.rows ();
00323   octave_idx_type b_nc = b.cols ();
00324   octave_idx_type b_nz = b.nnz ();
00325   octave_idx_type nz = 0;
00326   a = MSparse<RT> (b_nr, b_nc, b_nz);
00327   octave_sort<octave_idx_type> sort;
00328   octave_idx_type *ri = a.xridx();
00329   OCTAVE_LOCAL_BUFFER (RT, X, b_nr);
00330   a.xcidx(0) = 0;
00331   for (octave_idx_type j = 0; j < b_nc; j++)
00332     {
00333       for (octave_idx_type i = b.cidx(j); i < b.cidx(j+1); i++)
00334         {
00335           octave_quit ();
00336           octave_idx_type r = p [b.ridx (i)];
00337           X [r] = b.data (i);
00338           a.xridx(nz++) = p [b.ridx (i)];
00339         }
00340       sort.sort (ri + a.xcidx (j), nz - a.xcidx (j));
00341       for (octave_idx_type i = a.cidx (j); i < nz; i++)
00342         {
00343           octave_quit ();
00344           a.xdata (i) = X [a.xridx (i)];
00345         }
00346       a.xcidx(j+1) = nz;
00347     }
00348 }
00349 
00350 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00351 static void
00352 dmsolve_permute (MSparse<double> &a, const MSparse<double>& b,
00353                  const octave_idx_type *p);
00354 
00355 static void
00356 dmsolve_permute (MSparse<Complex> &a, const MSparse<double>& b,
00357                  const octave_idx_type *p);
00358 
00359 static void
00360 dmsolve_permute (MSparse<Complex> &a, const MSparse<Complex>& b,
00361                  const octave_idx_type *p);
00362 #endif
00363 
00364 static void
00365 solve_singularity_warning (double)
00366 {
00367   // Dummy singularity handler so that LU solver doesn't flag
00368   // an error for numerically rank defficient matrices
00369 }
00370 
00371 template <class RT, class ST, class T>
00372 RT
00373 dmsolve (const ST &a, const T &b, octave_idx_type &info)
00374 {
00375 #ifdef HAVE_CXSPARSE
00376   octave_idx_type nr = a.rows ();
00377   octave_idx_type nc = a.cols ();
00378   octave_idx_type b_nr = b.rows ();
00379   octave_idx_type b_nc = b.cols ();
00380   RT retval;
00381 
00382   if (nr < 0 || nc < 0 || nr != b_nr)
00383     (*current_liboctave_error_handler)
00384       ("matrix dimension mismatch in solution of minimum norm problem");
00385   else if (nr == 0 || nc == 0 || b_nc == 0)
00386     retval = RT (nc, b_nc, 0.0);
00387   else
00388     {
00389       octave_idx_type nnz_remaining = a.nnz ();
00390       CXSPARSE_DNAME () csm;
00391       csm.m = nr;
00392       csm.n = nc;
00393       csm.x = 0;
00394       csm.nz = -1;
00395       csm.nzmax = a.nnz ();
00396       // Cast away const on A, with full knowledge that CSparse won't touch it.
00397       // Prevents the methods below making a copy of the data.
00398       csm.p = const_cast<octave_idx_type *>(a.cidx ());
00399       csm.i = const_cast<octave_idx_type *>(a.ridx ());
00400 
00401 #if defined(CS_VER) && (CS_VER >= 2)
00402       CXSPARSE_DNAME (d) *dm = CXSPARSE_DNAME(_dmperm) (&csm, 0);
00403       octave_idx_type *p = dm->p;
00404       octave_idx_type *q = dm->q;
00405 #else
00406       CXSPARSE_DNAME (d) *dm = CXSPARSE_DNAME(_dmperm) (&csm);
00407       octave_idx_type *p = dm->P;
00408       octave_idx_type *q = dm->Q;
00409 #endif
00410       OCTAVE_LOCAL_BUFFER (octave_idx_type, pinv, nr);
00411       for (octave_idx_type i = 0; i < nr; i++)
00412         pinv [p [i]] = i;
00413       RT btmp;
00414       dmsolve_permute (btmp, b, pinv);
00415       info = 0;
00416       retval.resize (nc, b_nc);
00417 
00418       // Leading over-determined block
00419       if (dm->rr [2] < nr && dm->cc [3] < nc)
00420         {
00421           ST m = dmsolve_extract (a, pinv, q, dm->rr [2], nr, dm->cc [3], nc,
00422                                   nnz_remaining, true);
00423           nnz_remaining -= m.nnz();
00424           RT mtmp =
00425             qrsolve (m, dmsolve_extract (btmp, 0, 0, dm->rr[2], b_nr, 0,
00426                                          b_nc), info);
00427           dmsolve_insert (retval, mtmp, q, dm->cc [3], 0);
00428           if (dm->rr [2] > 0 && !info)
00429             {
00430               m = dmsolve_extract (a, pinv, q, 0, dm->rr [2],
00431                                    dm->cc [3], nc, nnz_remaining, true);
00432               nnz_remaining -= m.nnz();
00433               RT ctmp = dmsolve_extract (btmp, 0, 0, 0,
00434                                          dm->rr[2], 0, b_nc);
00435               btmp.insert (ctmp - m * mtmp, 0, 0);
00436             }
00437         }
00438 
00439       // Structurally non-singular blocks
00440       // FIXME Should use fine Dulmange-Mendelsohn decomposition here.
00441       if (dm->rr [1] < dm->rr [2] && dm->cc [2] < dm->cc [3] && !info)
00442         {
00443           ST m = dmsolve_extract (a, pinv, q, dm->rr [1], dm->rr [2],
00444                                   dm->cc [2], dm->cc [3], nnz_remaining, false);
00445           nnz_remaining -= m.nnz();
00446           RT btmp2 = dmsolve_extract (btmp, 0, 0, dm->rr [1], dm->rr [2],
00447                                       0, b_nc);
00448           double rcond = 0.0;
00449           MatrixType mtyp (MatrixType::Full);
00450           RT mtmp = m.solve (mtyp, btmp2, info, rcond,
00451                              solve_singularity_warning, false);
00452           if (info != 0)
00453             {
00454               info = 0;
00455               mtmp = qrsolve (m, btmp2, info);
00456             }
00457 
00458           dmsolve_insert (retval, mtmp, q, dm->cc [2], 0);
00459           if (dm->rr [1] > 0 && !info)
00460             {
00461               m = dmsolve_extract (a, pinv, q, 0, dm->rr [1], dm->cc [2],
00462                                    dm->cc [3], nnz_remaining, true);
00463               nnz_remaining -= m.nnz();
00464               RT ctmp = dmsolve_extract (btmp, 0, 0, 0,
00465                                          dm->rr[1], 0, b_nc);
00466               btmp.insert (ctmp - m * mtmp, 0, 0);
00467             }
00468         }
00469 
00470       // Trailing under-determined block
00471       if (dm->rr [1] > 0 && dm->cc [2] > 0 && !info)
00472         {
00473           ST m = dmsolve_extract (a, pinv, q, 0, dm->rr [1], 0,
00474                                   dm->cc [2], nnz_remaining, true);
00475           RT mtmp =
00476             qrsolve (m, dmsolve_extract(btmp, 0, 0, 0, dm->rr [1] , 0,
00477                                         b_nc), info);
00478           dmsolve_insert (retval, mtmp, q, 0, 0);
00479         }
00480 
00481       CXSPARSE_DNAME (_dfree) (dm);
00482     }
00483   return retval;
00484 #else
00485   return RT ();
00486 #endif
00487 }
00488 
00489 #if !defined (CXX_NEW_FRIEND_TEMPLATE_DECL)
00490 extern Matrix
00491 dmsolve (const SparseMatrix &a, const Matrix &b,
00492          octave_idx_type &info);
00493 
00494 extern ComplexMatrix
00495 dmsolve (const SparseMatrix &a, const ComplexMatrix &b,
00496          octave_idx_type &info);
00497 
00498 extern ComplexMatrix
00499 dmsolve (const SparseComplexMatrix &a, const Matrix &b,
00500          octave_idx_type &info);
00501 
00502 extern ComplexMatrix
00503 dmsolve (const SparseComplexMatrix &a, const ComplexMatrix &b,
00504          octave_idx_type &info);
00505 
00506 extern SparseMatrix
00507 dmsolve (const SparseMatrix &a, const SparseMatrix &b,
00508          octave_idx_type &info);
00509 
00510 extern SparseComplexMatrix
00511 dmsolve (const SparseMatrix &a, const SparseComplexMatrix &b,
00512          octave_idx_type &info);
00513 
00514 extern SparseComplexMatrix
00515 dmsolve (const SparseComplexMatrix &a, const SparseMatrix &b,
00516          octave_idx_type &info);
00517 
00518 extern SparseComplexMatrix
00519 dmsolve (const SparseComplexMatrix &a, const SparseComplexMatrix &b,
00520          octave_idx_type &info);
00521 #endif
 All Classes Files Functions Variables Typedefs Enumerations Enumerator Friends Defines