#include #include "intflapack.h" #include "tvector.h" #include "tmatrix.h" #include /*! \defgroup LinAlg LinAlg module This module contains classes and functions for complex linear algebra on arrays. This module is intended mainly to have classes implementing C++ interfaces between Sophya objects and external linear algebra libraries, such as LAPACK. */ /*! \class SOPHYA::LapackServer \ingroup LinAlg This class implements an interface to LAPACK library driver routines. The LAPACK (Linear Algebra PACKage) is a collection high performance routines to solve common problems in numerical linear algebra. its is available from http://www.netlib.org. The present version of our LapackServer (Feb 2001) provides only interfaces for the linear system solver and singular value decomposition (SVD). Only arrays with BaseArray::FortranMemoryMapping can be handled by LapackServer. LapackServer can be instanciated for simple and double precision real or complex array types. The example below shows solving a linear system A*X = B \code #include "intflapack.h" // ... // Use FortranMemoryMapping as default BaseArray::SetDefaultMemoryMapping(BaseArray::FortranMemoryMapping); // Create an fill the arrays A and B int n = 20; Matrix A(n, n); A = RandomSequence(); Vector X(n),B(n); X = RandomSequence(); B = A*X; // Solve the linear system A*X = B LapackServer lps; lps.LinSolve(A,B); // We get the result in B, which should be equal to X ... // Compute the difference B-X ; Vector diff = B-X; \endcode */ extern "C" { // Drivers pour resolution de systemes lineaires void sgesv_(int_4* n, int_4* nrhs, r_4* a, int_4* lda, int_4* ipiv, r_4* b, int_4* ldb, int_4* info); void dgesv_(int_4* n, int_4* nrhs, r_8* a, int_4* lda, int_4* ipiv, r_8* b, int_4* ldb, int_4* info); void cgesv_(int_4* n, int_4* nrhs, complex* a, int_4* lda, int_4* ipiv, complex* b, int_4* ldb, int_4* info); void zgesv_(int_4* n, int_4* nrhs, complex* a, int_4* lda, int_4* ipiv, complex* b, int_4* ldb, int_4* info); // Driver pour decomposition SVD void sgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, r_4* a, int_4* lda, r_4* s, r_4* u, int_4* ldu, r_4* vt, int_4* ldvt, r_4* work, int_4* lwork, int_4* info); void dgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, r_8* a, int_4* lda, r_8* s, r_8* u, int_4* ldu, r_8* vt, int_4* ldvt, r_8* work, int_4* lwork, int_4* info); void cgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, complex* a, int_4* lda, complex* s, complex* u, int_4* ldu, complex* vt, int_4* ldvt, complex* work, int_4* lwork, int_4* info); void zgesvd_(char* jobu, char* jobvt, int_4* m, int_4* n, complex* a, int_4* lda, complex* s, complex* u, int_4* ldu, complex* vt, int_4* ldvt, complex* work, int_4* lwork, int_4* info); } // -------------- Classe LapackServer -------------- template LapackServer::LapackServer() { SetWorkSpaceSizeFactor(); } template LapackServer::~LapackServer() { } //! Interface to Lapack linear system solver driver s/d/c/zgesvd(). /*! Solve the linear system a * x = b. Input arrays should have FortranMemory mapping (column packed). \param a : input matrix, overwritten on output \param b : input-output, input vector b, contains x on exit \return : return code from lapack driver _gesv() */ template int LapackServer::LinSolve(TArray& a, TArray & b) { if ( ( a.NbDimensions() != 2 ) || ( b.NbDimensions() != 2 ) ) throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b NbDimensions() != 2")); int_4 rowa = a.RowsKA(); int_4 cola = a.ColsKA(); int_4 rowb = b.RowsKA(); int_4 colb = b.ColsKA(); if ( a.Size(rowa) != a.Size(cola)) throw(SzMismatchError("LapackServer::LinSolve(a,b) a Not a square Array")); if ( a.Size(rowa) != b.Size(rowb)) throw(SzMismatchError("LapackServer::LinSolve(a,b) RowSize(a <> b) ")); if (!a.IsPacked(rowa) || !b.IsPacked(rowb)) throw(SzMismatchError("LapackServer::LinSolve(a,b) a Or b Not Column Packed")); int_4 n = a.Size(rowa); int_4 nrhs = b.Size(colb); int_4 lda = a.Step(cola); int_4 ldb = b.Step(colb); int_4 info; int_4* ipiv = new int_4[n]; if (typeid(T) == typeid(r_4) ) sgesv_(&n, &nrhs, (r_4 *)a.Data(), &lda, ipiv, (r_4 *)b.Data(), &ldb, &info); else if (typeid(T) == typeid(r_8) ) dgesv_(&n, &nrhs, (r_8 *)a.Data(), &lda, ipiv, (r_8 *)b.Data(), &ldb, &info); else if (typeid(T) == typeid(complex) ) cgesv_(&n, &nrhs, (complex *)a.Data(), &lda, ipiv, (complex *)b.Data(), &ldb, &info); else if (typeid(T) == typeid(complex) ) zgesv_(&n, &nrhs, (complex *)a.Data(), &lda, ipiv, (complex *)b.Data(), &ldb, &info); else { delete[] ipiv; string tn = typeid(T).name(); cerr << " LapackServer::LinSolve(a,b) - Unsupported DataType T = " << tn << endl; throw TypeMismatchExc("LapackServer::LinSolve(a,b) - Unsupported DataType (T)"); } delete[] ipiv; return(info); } //! Interface to Lapack SVD driver s/d/c/zgesv(). /*! Computes the vector of singular values of \b a. Input arrays should have FortranMemoryMapping (column packed). \param a : input m-by-n matrix \param s : Vector of min(m,n) singular values (descending order) \return : return code from lapack driver _gesvd() */ template int LapackServer::SVD(TArray& a, TArray & s) { return (SVDDriver(a, s, NULL, NULL) ); } //! Interface to Lapack SVD driver s/d/c/zgesv(). /*! Computes the vector of singular values of \b a, as well as right and left singular vectors of \b a. \f[ A = U \Sigma V^T , ( A = U \Sigma V^H \ complex) \f] \f[ A v_i = \sigma_i u_i \ and A^T u_i = \sigma_i v_i \ (A^H \ complex) \f] U and V are orthogonal (unitary) matrices. \param a : input m-by-n matrix (in FotranMemoryMapping) \param s : Vector of min(m,n) singular values (descending order) \param u : Matrix of left singular vectors \param vt : Transpose of right singular vectors. \return : return code from lapack driver _gesvd() */ template int LapackServer::SVD(TArray& a, TArray & s, TArray & u, TArray & vt) { return (SVDDriver(a, s, &u, &vt) ); } //! Interface to Lapack SVD driver s/d/c/zgesv(). template int LapackServer::SVDDriver(TArray& a, TArray & s, TArray* up, TArray* vtp) { if ( ( a.NbDimensions() != 2 ) ) throw(SzMismatchError("LapackServer::SVD(a, ...) a.NbDimensions() != 2")); int_4 rowa = a.RowsKA(); int_4 cola = a.ColsKA(); if ( !a.IsPacked(rowa) ) throw(SzMismatchError("LapackServer::SVD(a, ...) a Not Column Packed ")); int_4 m = a.Size(rowa); int_4 n = a.Size(cola); int_4 maxmn = (m > n) ? m : n; int_4 minmn = (m < n) ? m : n; char jobu, jobvt; jobu = 'N'; jobvt = 'N'; sa_size_t sz[2]; if ( up != NULL) { if ( dynamic_cast< TVector * > (vtp) ) throw( TypeMismatchExc("LapackServer::SVD() Wrong type (=TVector) for u !") ); up->SetMemoryMapping(BaseArray::FortranMemoryMapping); sz[0] = sz[1] = m; up->ReSize(2, sz ); jobu = 'A'; } else { up = new TMatrix(1,1); jobu = 'N'; } if ( vtp != NULL) { if ( dynamic_cast< TVector * > (vtp) ) throw( TypeMismatchExc("LapackServer::SVD() Wrong type (=TVector) for vt !") ); vtp->SetMemoryMapping(BaseArray::FortranMemoryMapping); sz[0] = sz[1] = n; vtp->ReSize(2, sz ); jobvt = 'A'; } else { vtp = new TMatrix(1,1); jobvt = 'N'; } TVector *vs = dynamic_cast< TVector * > (&s); if (vs) vs->ReSize(minmn); else { TMatrix *ms = dynamic_cast< TMatrix * > (&s); if (ms) ms->ReSize(minmn,1); else { sz[0] = minmn; sz[1] = 1; s.ReSize(1, sz); } } int_4 lda = a.Step(a.ColsKA()); int_4 ldu = up->Step(up->ColsKA()); int_4 ldvt = vtp->Step(vtp->ColsKA()); int_4 lwork = maxmn*5*wspace_size_factor; T * work = new T[lwork]; int_4 info; if (typeid(T) == typeid(r_4) ) sgesvd_(&jobu, &jobvt, &m, &n, (r_4 *)a.Data(), &lda, (r_4 *)s.Data(), (r_4 *) up->Data(), &ldu, (r_4 *)vtp->Data(), &ldvt, (r_4 *)work, &lwork, &info); else if (typeid(T) == typeid(r_8) ) dgesvd_(&jobu, &jobvt, &m, &n, (r_8 *)a.Data(), &lda, (r_8 *)s.Data(), (r_8 *) up->Data(), &ldu, (r_8 *)vtp->Data(), &ldvt, (r_8 *)work, &lwork, &info); else if (typeid(T) == typeid(complex) ) cgesvd_(&jobu, &jobvt, &m, &n, (complex *)a.Data(), &lda, (complex *)s.Data(), (complex *) up->Data(), &ldu, (complex *)vtp->Data(), &ldvt, (complex *)work, &lwork, &info); else if (typeid(T) == typeid(complex) ) zgesvd_(&jobu, &jobvt, &m, &n, (complex *)a.Data(), &lda, (complex *)s.Data(), (complex *) up->Data(), &ldu, (complex *)vtp->Data(), &ldvt, (complex *)work, &lwork, &info); else { if (jobu == 'N') delete up; if (jobvt == 'N') delete vtp; string tn = typeid(T).name(); cerr << " LapackServer::SVDDriver(...) - Unsupported DataType T = " << tn << endl; throw TypeMismatchExc("LapackServer::LinSolve(a,b) - Unsupported DataType (T)"); } if (jobu == 'N') delete up; if (jobvt == 'N') delete vtp; return(info); } void rztest_lapack(TArray& aa, TArray& bb) { if ( aa.NbDimensions() != 2 ) throw(SzMismatchError("rztest_lapack(TMatrix A Not a Matrix")); if ( aa.SizeX() != aa.SizeY()) throw(SzMismatchError("rztest_lapack(TMatrix A Not a square Matrix")); if ( bb.NbDimensions() != 2 ) throw(SzMismatchError("rztest_lapack(TMatrix A Not a Matrix")); if ( bb.SizeX() != aa.SizeX() ) throw(SzMismatchError("rztest_lapack(TMatrix A <> B ")); if ( !bb.IsPacked() || !bb.IsPacked() ) throw(SzMismatchError("rztest_lapack(TMatrix Not packed A or B ")); int_4 n = aa.SizeX(); int_4 nrhs = bb.SizeY(); int_4 lda = n; int_4 ldb = bb.SizeX(); int_4 info; int_4* ipiv = new int_4[n]; sgesv_(&n, &nrhs, aa.Data(), &lda, ipiv, bb.Data(), &ldb, &info); delete[] ipiv; cout << "rztest_lapack/Info= " << info << endl; cout << aa << "\n" << bb << endl; return; } /////////////////////////////////////////////////////////////// #ifdef __CXX_PRAGMA_TEMPLATES__ #pragma define_template LapackServer #pragma define_template LapackServer #pragma define_template LapackServer< complex > #pragma define_template LapackServer< complex > #endif #if defined(ANSI_TEMPLATES) || defined(GNU_TEMPLATES) template class LapackServer; template class LapackServer; template class LapackServer< complex >; template class LapackServer< complex >; #endif #if defined(OS_LINUX) // Pour le link avec f2c sous Linux extern "C" { void MAIN__(); } void MAIN__() { cerr << "MAIN__() function for linking with libf2c.a " << endl; cerr << " This function should never be called !!! " << endl; throw PError("MAIN__() should not be called - see intflapack.cc"); } #endif