SNAP Library 2.2, Developer Reference  2014-03-11 19:15:55
SNAP, a general purpose, high performance system for analysis and manipulation of large networks
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
linalg.h
Go to the documentation of this file.
00001 
00002 // forward declarations
00003 class TLinAlg;
00004 class TLAMisc;
00005 
00007 // Matrix
00008 class TMatrix {
00009 private:
00010     bool Transposed;
00011 protected:
00012     virtual void PMultiply(const TFltVV& B, int ColId, TFltV& Result) const = 0;
00013     virtual void PMultiply(const TFltV& Vec, TFltV& Result) const = 0;
00014     virtual void PMultiplyT(const TFltVV& B, int ColId, TFltV& Result) const = 0;
00015     virtual void PMultiplyT(const TFltV& Vec, TFltV& Result) const = 0;
00016 
00017     virtual int PGetRows() const = 0;
00018     virtual int PGetCols() const = 0;
00019 public:
00020     TMatrix(): Transposed(false) {}
00021     virtual ~TMatrix() { }
00022 
00023     // Result = A * B(:,ColId)
00024     void Multiply(const TFltVV& B, int ColId, TFltV& Result) const {
00025         if (Transposed) { PMultiplyT(B, ColId, Result); }
00026         else { PMultiply(B, ColId, Result); }
00027     }
00028     // Result = A * Vec
00029     void Multiply(const TFltV& Vec, TFltV& Result) const {
00030         if (Transposed) { PMultiplyT(Vec, Result); }
00031         else { PMultiply(Vec, Result); }
00032     }
00033     // Result = A' * B(:,ColId)
00034     void MultiplyT(const TFltVV& B, int ColId, TFltV& Result) const {
00035         if (Transposed) { PMultiply(B, ColId, Result); }
00036         else { PMultiplyT(B, ColId, Result); }
00037     }
00038     // Result = A' * Vec
00039     void MultiplyT(const TFltV& Vec, TFltV& Result) const{
00040         if (Transposed) { PMultiply(Vec, Result); }
00041         else { PMultiplyT(Vec, Result); }
00042     }
00043 
00044     // number of rows
00045     int GetRows() const { return Transposed ? PGetCols() : PGetRows(); }
00046     // number of columns
00047     int GetCols() const { return Transposed ? PGetRows() : PGetCols(); }
00048 
00049     void Transpose() { Transposed = !Transposed; }
00050 };
00051 
00053 // Sparse-Column-Matrix
00054 //  matrix is given with columns as sparse vectors
00055 class TSparseColMatrix: public TMatrix {
00056 public:
00057     // number of rows and columns of matrix
00058     int RowN, ColN;
00059     // vector of sparse columns
00060     TVec<TIntFltKdV> ColSpVV;
00061 protected:
00062     // Result = A * B(:,ColId)
00063     virtual void PMultiply(const TFltVV& B, int ColId, TFltV& Result) const;
00064     // Result = A * Vec
00065     virtual void PMultiply(const TFltV& Vec, TFltV& Result) const;
00066     // Result = A' * B(:,ColId)
00067     virtual void PMultiplyT(const TFltVV& B, int ColId, TFltV& Result) const;
00068     // Result = A' * Vec
00069     virtual void PMultiplyT(const TFltV& Vec, TFltV& Result) const;
00070 
00071     int PGetRows() const { return RowN; }
00072     int PGetCols() const { return ColN; }
00073 
00074 public:
00075     TSparseColMatrix(): TMatrix() {}
00076     TSparseColMatrix(TVec<TIntFltKdV> _ColSpVV): TMatrix(), ColSpVV(_ColSpVV) {}
00077     TSparseColMatrix(TVec<TIntFltKdV> _ColSpVV, const int& _RowN, const int& _ColN): 
00078                 TMatrix(), RowN(_RowN), ColN(_ColN), ColSpVV(_ColSpVV) {}
00079     // loads Matlab sparse matrix format: row, column, value.
00080     //   Indexes start with 1.
00081     void Save(TSOut& SOut) {
00082         SOut.Save(RowN); SOut.Save(ColN); ColSpVV.Save(SOut); }
00083     void Load(TSIn& SIn) {
00084         SIn.Load(RowN); SIn.Load(ColN); ColSpVV = TVec<TIntFltKdV>(SIn); }
00085 };
00086 
00088 // Sparse-Row-Matrix
00089 //  matrix is given with rows as sparse vectors
00090 class TSparseRowMatrix: public TMatrix {
00091 public:
00092     // number of rows and columns of matrix
00093     int RowN, ColN;
00094     // vector of sparse rows
00095     TVec<TIntFltKdV> RowSpVV;
00096 protected:
00097     // Result = A * B(:,ColId)
00098     virtual void PMultiply(const TFltVV& B, int ColId, TFltV& Result) const;
00099     // Result = A * Vec
00100     virtual void PMultiply(const TFltV& Vec, TFltV& Result) const;
00101     // Result = A' * B(:,ColId)
00102     virtual void PMultiplyT(const TFltVV& B, int ColId, TFltV& Result) const;
00103     // Result = A' * Vec
00104     virtual void PMultiplyT(const TFltV& Vec, TFltV& Result) const;
00105 
00106     int PGetRows() const { return RowN; }
00107     int PGetCols() const { return ColN; }
00108 
00109 public:
00110     TSparseRowMatrix(): TMatrix() {}
00111     TSparseRowMatrix(TVec<TIntFltKdV> _RowSpVV): TMatrix(), RowSpVV(_RowSpVV) {}
00112     TSparseRowMatrix(TVec<TIntFltKdV> _RowSpVV, const int& _RowN, const int& _ColN): 
00113                 TMatrix(), RowN(_RowN), ColN(_ColN), RowSpVV(_RowSpVV) {}
00114         // loads Matlab sparse matrix format: row, column, value.
00115     //   Indexes start with 1.
00116     TSparseRowMatrix(const TStr& MatlabMatrixFNm);
00117     void Save(TSOut& SOut) {
00118         SOut.Save(RowN); SOut.Save(ColN); RowSpVV.Save(SOut); }
00119     void Load(TSIn& SIn) {
00120         SIn.Load(RowN); SIn.Load(ColN); RowSpVV = TVec<TIntFltKdV>(SIn); }
00121 };
00122 
00124 // Full-Col-Matrix
00125 //  matrix is given with columns of full vectors
00126 class TFullColMatrix: public TMatrix {
00127 public:
00128     // number of rows and columns of matrix
00129     int RowN, ColN;
00130     // vector of sparse columns
00131     TVec<TFltV> ColV;
00132 protected:
00133     // Result = A * B(:,ColId)
00134     virtual void PMultiply(const TFltVV& B, int ColId, TFltV& Result) const;
00135     // Result = A * Vec
00136     virtual void PMultiply(const TFltV& Vec, TFltV& Result) const;
00137     // Result = A' * B(:,ColId)
00138     virtual void PMultiplyT(const TFltVV& B, int ColId, TFltV& Result) const;
00139     // Result = A' * Vec
00140     virtual void PMultiplyT(const TFltV& Vec, TFltV& Result) const;
00141 
00142     int PGetRows() const { return RowN; }
00143     int PGetCols() const { return ColN; }
00144 
00145 public:
00146     TFullColMatrix(): TMatrix() {}
00147     // loads matrix saved in matlab with command:
00148     //  save -ascii Matrix.dat M
00149     TFullColMatrix(const TStr& MatlabMatrixFNm);
00150     void Save(TSOut& SOut) { ColV.Save(SOut); }
00151     void Load(TSIn& SIn) { ColV.Load(SIn); }
00152 };
00153 
00155 // Basic Linear Algebra Operations
00156 class TLinAlg {
00157 public:
00158     // <x,y>
00159     static double DotProduct(const TFltV& x, const TFltV& y);
00160     // <X(:,ColIdX), Y(:,ColIdY)>
00161     static double DotProduct(const TFltVV& X, int ColIdX, const TFltVV& Y, int ColIdY);
00162     // <X(:,ColId), Vec>
00163     static double DotProduct(const TFltVV& X, int ColId, const TFltV& Vec);
00164     // sparse dot products:
00165     // <x,y> where x AND y are sparse
00166     static double DotProduct(const TIntFltKdV& x, const TIntFltKdV& y);
00167     // <x,y> where only y is sparse
00168     static double DotProduct(const TFltV& x, const TIntFltKdV& y);
00169     // <X(:,ColId),y> where only y is sparse
00170     static double DotProduct(const TFltVV& X, int ColId, const TIntFltKdV& y);
00171 
00172     // z := p * x + q * y
00173     static void LinComb(const double& p, const TFltV& x,
00174         const double& q, const TFltV& y, TFltV& z);
00175     // z := p * x + (1 - p) * y
00176     static void ConvexComb(const double& p, const TFltV& x, const TFltV& y, TFltV& z);
00177 
00178     // z := k * x + y
00179     static void AddVec(const double& k, const TFltV& x, const TFltV& y, TFltV& z);
00180     // z := k * x + y
00181     static void AddVec(const double& k, const TIntFltKdV& x, const TFltV& y, TFltV& z);
00182     // y := k * x + y
00183     static void AddVec(const double& k, const TIntFltKdV& x, TFltV& y);
00184     // Y(:,Col) += k * X(:,Col)
00185     static void AddVec(double k, const TFltVV& X, int ColIdX, TFltVV& Y, int ColIdY);
00186     // Result += k * X(:,Col)
00187     static void AddVec(double k, const TFltVV& X, int ColId, TFltV& Result);
00188         // z = x + y
00189     static void AddVec(const TIntFltKdV& x, const TIntFltKdV& y, TIntFltKdV& z);
00190 
00191     // Result = SUM(x)
00192     static double SumVec(const TFltV& x);
00193     // Result = SUM(k*x + y)
00194     static double SumVec(double k, const TFltV& x, const TFltV& y);
00195 
00196     // Result = ||x-y||^2 (Euclidian)
00197     static double EuclDist2(const TFltV& x, const TFltV& y);
00198     static double EuclDist2(const TFltPr& x, const TFltPr& y);
00199     // Result = ||x-y|| (Euclidian)
00200     static double EuclDist(const TFltV& x, const TFltV& y);
00201     static double EuclDist(const TFltPr& x, const TFltPr& y);
00202 
00203     // ||x||^2 (Euclidian)
00204     static double Norm2(const TFltV& x);
00205     // ||x|| (Euclidian)
00206     static double Norm(const TFltV& x);
00207     // x := x / ||x||
00208     static void Normalize(TFltV& x);
00209 
00210     // ||x||^2 (Euclidian), x is sparse
00211     static double Norm2(const TIntFltKdV& x);
00212     // ||x|| (Euclidian), x is sparse
00213     static double Norm(const TIntFltKdV& x);
00214     // x := x / ||x||, x is sparse
00215     static void Normalize(TIntFltKdV& x);
00216 
00217     // ||X(:,ColId)||^2 (Euclidian)
00218     static double Norm2(const TFltVV& X, int ColId);
00219     // ||X(:,ColId)|| (Euclidian)
00220     static double Norm(const TFltVV& X, int ColId);
00221 
00222     // L1 norm of x (Sum[|xi|, i = 1..n])
00223     static double NormL1(const TFltV& x);
00224     // L1 norm of k*x+y (Sum[|k*xi+yi|, i = 1..n])
00225     static double NormL1(double k, const TFltV& x, const TFltV& y);
00226     // L1 norm of x (Sum[|xi|, i = 1..n])
00227     static double NormL1(const TIntFltKdV& x);
00228     // x := x / ||x||_inf
00229     static void NormalizeL1(TFltV& x);
00230     // x := x / ||x||_inf
00231     static void NormalizeL1(TIntFltKdV& x);
00232 
00233     // Linf norm of x (Max{|xi|, i = 1..n})
00234     static double NormLinf(const TFltV& x);
00235     // Linf norm of x (Max{|xi|, i = 1..n})
00236     static double NormLinf(const TIntFltKdV& x);
00237     // x := x / ||x||_inf
00238     static void NormalizeLinf(TFltV& x);
00239     // x := x / ||x||_inf, , x is sparse
00240     static void NormalizeLinf(TIntFltKdV& x);
00241 
00242     // y := k * x
00243     static void MultiplyScalar(const double& k, const TFltV& x, TFltV& y);
00244     // y := k * x
00245     static void MultiplyScalar(const double& k, const TIntFltKdV& x, TIntFltKdV& y);
00246 
00247     // y := A * x
00248     static void Multiply(const TFltVV& A, const TFltV& x, TFltV& y);
00249     // C(:, ColId) := A * x
00250     static void Multiply(const TFltVV& A, const TFltV& x, TFltVV& C, int ColId);
00251     // y := A * B(:, ColId)
00252     static void Multiply(const TFltVV& A, const TFltVV& B, int ColId, TFltV& y);
00253     // C(:, ColIdC) := A * B(:, ColIdB)
00254     static void Multiply(const TFltVV& A, const TFltVV& B, int ColIdB, TFltVV& C, int ColIdC);
00255 
00256     // y := A' * x
00257     static void MultiplyT(const TFltVV& A, const TFltV& x, TFltV& y);
00258 
00259     // C = A * B
00260     static void Multiply(const TFltVV& A, const TFltVV& B, TFltVV& C);
00261         
00262         // D = alpha * A(') * B(') + beta * C(')
00263         typedef enum { GEMM_NO_T = 0, GEMM_A_T = 1, GEMM_B_T = 2, GEMM_C_T = 4 } TLinAlgGemmTranspose;
00264         static void Gemm(const double& Alpha, const TFltVV& A, const TFltVV& B, const double& Beta, 
00265                 const TFltVV& C, TFltVV& D, const int& TransposeFlags);
00266         
00267         // B = A^(-1)
00268         typedef enum { DECOMP_SVD } TLinAlgInverseType;
00269         static void Inverse(const TFltVV& A, TFltVV& B, const TLinAlgInverseType& DecompType);
00270         // subtypes of finding an inverse
00271         static void InverseSVD(const TFltVV& A, TFltVV& B);
00272 
00273         // transpose matrix - B = A'
00274         static void Transpose(const TFltVV& A, TFltVV& B);
00275 
00276     // performes Gram-Schmidt ortogonalization on elements of Q
00277     static void GS(TVec<TFltV>& Q);
00278     // Gram-Schmidt on columns of matrix Q
00279     static void GS(TFltVV& Q);
00280 
00281     // rotates vector (OldX,OldY) for angle Angle (in radians!)
00282     static void Rotate(const double& OldX, const double& OldY, const double& Angle, double& NewX, double& NewY);
00283 
00284     // checks if set of vectors is ortogonal
00285     static void AssertOrtogonality(const TVec<TFltV>& Vecs, const double& Threshold);
00286     static void AssertOrtogonality(const TFltVV& Vecs, const double& Threshold);
00287 };
00288 
00290 // Numerical-Recepies-Exception
00291 class TNSException : public TExcept {
00292 public:
00293     TStr Message;
00294 public:
00295     TNSException(const TStr& Msg): TExcept(Msg) {}
00296 };
00297 
00299 // Numerical-Linear-Algebra (copied from Numerical Recepies)
00300 class TNumericalStuff {
00301 private:
00302   static double sqr(double a);
00303   static double sign(double a, double b);
00304 
00305   // Computes (a^2 + b^2)^(1/2) without
00306   // destructive underflow or overflow.
00307   static double pythag(double a, double b);
00308 
00309   //displays error message to screen
00310   static void nrerror(const TStr& error_text);
00311 
00312 public:
00313     // Householder reduction of a real, symmetric matrix a[1..n][1..n].
00314     // On output, a is replaced by the orthogonal matrix Q eecting the
00315     // transformation. d[1..n] returns the diagonal elements of the
00316     // tridiagonal matrix, and e[1..n] the o-diagonal elements, with
00317     // e[1]=0. Several statements, as noted in comments, can be omitted
00318     // if only eigenvalues are to be found, in which case a contains no
00319     // useful information on output. Otherwise they are to be included.
00320     static void SymetricToTridiag(TFltVV& a, int n, TFltV& d, TFltV& e);
00321 
00322         // QL algorithm with implicit shifts, to determine the eigenvalues
00323         // and eigenvectors of a real, symmetric, tridiagonal matrix, or of
00324         // a real, symmetric matrix previously reduced by tred2 x11.2. On
00325         // input, d[1..n] contains the diagonal elements of the tridiagonal
00326         // matrix. On output, it returns the eigenvalues. The vector e[1..n]
00327         // inputs the subdiagonal elements of the tridiagonal matrix, with
00328         // e[1] arbitrary. On output e is destroyed. When finding only the
00329         // eigenvalues, several lines may be omitted, as noted in the comments.
00330         // If the eigenvectors of a tridiagonal matrix are desired, the matrix
00331         // z[1..n][1..n] is input as the identity matrix. If the eigenvectors
00332         // of a matrix that has been reduced by tred2 are required, then z is
00333         // input as the matrix output by tred2. In either case, the kth column
00334         // of z returns the normalized eigenvector corresponding to d[k].
00335         static void EigSymmetricTridiag(TFltV& d, TFltV& e, int n, TFltVV& z);
00336 
00337         // Given a positive-dedinite symmetric matrix A(n,n), this routine
00338         // constructs its Cholesky decomposition, A = L * L^T . On input, only
00339         // the upper triangle of A need be given; it is not modified. The
00340         // Cholesky factor L is returned in the lower triangle of A, except for
00341         // its diagonal elements which are returned in p(n).
00342         static void CholeskyDecomposition(TFltVV& A, TFltV& p);
00343 
00344         // Solves the set of n linear equations A * x = b, where A is a
00345         // positive-definite symmetric matrix. A(n,n) and p[1..n] are input
00346         // as the output of the routine choldc. Only the lower triangle of A
00347         // is accessed. b(n) is input as the right-hand side vector. The
00348         // solution vector is returned in x(n). A  and p are not modified and
00349         // can be left in place for successive calls with diferent right-hand
00350         // sides b. b is not modified unless you identify b and x in the calling
00351         // sequence, which is allowed.
00352         static void CholeskySolve(const TFltVV& A, const TFltV& p, const TFltV& b, TFltV& x);
00353 
00354         // Solves system of linear equations A * x = b, where A is symetric
00355         // positive-definite matrix. A is first decomposed using
00356         // CholeskyDecomposition and after solved using CholeskySolve. Only
00357         // upper triangle of A need be given and it is not modified. However,
00358         // lower triangle is modified!
00359         static void SolveSymetricSystem(TFltVV& A, const TFltV& b, TFltV& x);
00360 
00361     // solve system A x_i = e_i for i = 1..n, where A and p are output
00362     // from CholeskyDecomposition. Result is stored to upper triangule
00363     // (possible since inverse of symetric matrix is also symetric! Sigh...)
00364     static void InverseSubstitute(TFltVV& A, const TFltV& p);
00365 
00366     // Calculates inverse of symetric positiv definit matrix
00367     // Matrix is given as upper triangule of A, result is stored
00368     // in upper triangule of A. Lower triangule is random (actually
00369     // it has part of Choleksy decompositon of A)
00370     static void InverseSymetric(TFltVV& A);
00371 
00372     // calcualtes inverse of upper triagonal matrix A
00373     // lower triangle is messed up...
00374     static void InverseTriagonal(TFltVV& A);
00375 
00376     // Given a matrix a[1..n][1..n], this routine replaces it by the LU
00377     // decomposition of a rowwise permutation of itself. a and n are input.
00378     // a is output, arranged as in equation (2.3.14) above; indx[1..n] is
00379     // an output vector that records the row permutation efected by the partial
00380     // pivoting; d is output as +-1 depending on whether the number of row
00381     // interchanges was even or odd, respectively. This routine is used in
00382     // combination with lubksb to solve linear equations or invert a matrix.
00383     static void LUDecomposition(TFltVV& A, TIntV& indx, double& d);
00384 
00385     // Solves the set of n linear equations A*X = B. Here a[1..n][1..n] is input,
00386     // not as the matrix A but rather as its LU decomposition, determined by the
00387     // routine ludcmp. indx[1..n] is input as the permutation vector returned by
00388     // ludcmp. b[1..n] is input as the right-hand side vector B, and returns with
00389     // the solution vector X. a, n, and indx are not modified by this routine and
00390     // can be left in place for successive calls with diferent right-hand sides b.
00391     // This routine takes into account the possibility that b will begin with many
00392     // zero elements, so it is efficient for use in matrix inversion.
00393     static void LUSolve(const TFltVV& A, const TIntV& indx, TFltV& b);
00394 
00395     // Solves system of linear equations A * x = b. A is first decomposed using
00396     // LUDecomposition and after solved using LUSolve. A is modified!
00397     static void SolveLinearSystem(TFltVV& A, const TFltV& b, TFltV& x);
00398 };
00399 
00401 // Sparse-SVD
00402 //   Calculates singular-value-decompositon for sparse matrixes.
00403 //   If A is a matrix than A is decomposed to A = U S V'
00404 //   where S is diagonal with singular values on diagonal and U
00405 //   and V are ortogonal (U'*U = V'*V = I).
00406 typedef enum { ssotNoOrto, ssotSelective, ssotFull } TSpSVDReOrtoType;
00407 class TSparseSVD {
00408 private:
00409     // Result = Matrix' * Matrix * Vec(:,ColId)
00410     static void MultiplyATA(const TMatrix& Matrix,
00411         const TFltVV& Vec, int ColId, TFltV& Result);
00412     // Result = Matrix' * Matrix * Vec
00413     static void MultiplyATA(const TMatrix& Matrix,
00414         const TFltV& Vec, TFltV& Result);
00415 public:
00416     // calculates NumEig eigen values of symetric matrix
00417     // if SvdMatrixProductP than matrix Matrix'*Matrix is used
00418     static void SimpleLanczos(const TMatrix& Matrix,
00419         const int& NumEig, TFltV& EigValV,
00420         const bool& DoLocalReortoP = false,
00421         const bool& SvdMatrixProductP = false);
00422     // fast, calculates NumEig largers eigen values and vectors
00423     // kk should be something like 4*NumEig
00424     // if SvdMatrixProductP than matrix Matrix'*Matrix is used
00425     static void Lanczos(const TMatrix& Matrix,
00426         int NumEig, int Iters, const TSpSVDReOrtoType& ReOrtoType,
00427         TFltV& EigValV, TFltVV& EigVecVV,
00428         const bool& SvdMatrixProductP = false);
00429     static void Lanczos2(const TMatrix& Matrix,
00430         int MaxNumEig, int MaxSecs, const TSpSVDReOrtoType& ReOrtoType,
00431         TFltV& EigValV, TFltVV& EigVecVV,
00432         const bool& SvdMatrixProductP = false);
00433 
00434     // calculates only singular values (based on SimpleLanczos)
00435     static void SimpleLanczosSVD(const TMatrix& Matrix,
00436         const int& CalcSV, TFltV& SngValV,
00437         const bool& DoLocalReortoP = false);
00438     // fast, calculates NumSV largers SV (based on Lanczos)
00439     static void LanczosSVD(const TMatrix& Matrix,
00440         int NumSV, int Iters, const TSpSVDReOrtoType& ReOrtoType,
00441         TFltV& SgnValV, TFltVV& LeftSgnVecVV, TFltVV& RightSgnVecVV);
00442 
00443     // slow - ortogonal iteration
00444     static void OrtoIterSVD(const TMatrix& Matrix,
00445         int NumSV, int IterN, TFltV& SgnValV);
00446 
00447     // projects sparse vector to space spanned by columns of matrix U
00448     static void Project(const TIntFltKdV& Vec, const TFltVV& U, TFltV& ProjVec);
00449 };
00450 
00452 // Sigmoid  --  made by Janez(TM)
00453 //  (y = 1/[1 + exp[-Ax+B]])
00454 class TSigmoid {
00455 private:
00456     TFlt A;
00457     TFlt B;
00458 private:
00459   // Evaluates how well the sigmoid function fits the data.
00460   // J(A, B) = - ln prod_i P(Y = y_i | Z = z_i).  The 'data' parameter
00461   // should contain (z_i, y_i) pairs.  Smaller J means a better fit.
00462   static double EvaluateFit(const TFltIntKdV& data, const double A, const double B);
00463   // Computes not only J but also its partial derivatives.
00464   static void EvaluateFit(const TFltIntKdV& data, const double A,
00465         const double B, double& J, double& JA, double& JB);
00466   // Let J(lambda) = J(A + lambda U, B + lambda V).
00467     // This function computes J and its first and second derivatives.
00468   // They can be used to choose a good lambda (using Newton's method)
00469     // when minimizing J. -- This method has not been tested yet.
00470   static void EvaluateFit(const TFltIntKdV& data, const double A,
00471         const double B, const double U, const double V, const double lambda,
00472     double& J, double& JJ, double& JJJ);
00473 public:
00474     TSigmoid() { };
00475     TSigmoid(const double& A_, const double& B_): A(A_), B(B_) { };
00476         // Tries to find a pair (A, B) that minimizes J(A, B).
00477     // Uses gradient descent.
00478     TSigmoid(const TFltIntKdV& data);
00479 
00480     TSigmoid(TSIn& SIn) { A.Load(SIn); B.Load(SIn); }
00481     void Load(TSIn& SIn) { A.Load(SIn); B.Load(SIn); }
00482     void Save(TSOut& SOut) const {A.Save(SOut); B.Save(SOut);}
00483 
00484     double GetVal(const double& x) const {
00485         return 1.0 / (1.0 + exp(-A * x + B)); }
00486     double operator()(const double& x) const {
00487         return GetVal(x); }
00488 
00489     void GetSigmoidAB(double& A_, double& B_) { A_=A; B_=B; }
00490 };
00491 
00493 // Useful stuff (hopefuly)
00494 class TLAMisc {
00495 public:
00496         // Dumps vector to file so Excel can read it
00497     static void SaveCsvTFltV(const TFltV& Vec, TSOut& SOut);
00498         // Dumps sparse vector to file so Matlab can read it
00499     static void SaveMatlabTFltIntKdV(const TIntFltKdV& SpV, const int& ColN, TSOut& SOut);
00500         // Dumps vector to file so Matlab can read it
00501     static void SaveMatlabTFltV(const TFltV& m, const TStr& FName);
00502         // Dumps vector to file so Matlab can read it
00503     static void SaveMatlabTIntV(const TIntV& m, const TStr& FName);
00504         // Dumps column ColId from m to file so Matlab can read it
00505     static void SaveMatlabTFltVVCol(const TFltVV& m, int ColId, const TStr& FName);
00506         // Dumps matrix to file so Matlab can read it
00507     static void SaveMatlabTFltVV(const TFltVV& m, const TStr& FName);
00508         // Dumps main minor rowN x colN to file so Matlab can read it
00509         static void SaveMatlabTFltVVMjrSubMtrx(const TFltVV& m, int rowN, int colN, const TStr& FName);
00510     // loads matlab full matrix
00511     static void LoadMatlabTFltVV(const TStr& FNm, TVec<TFltV>& ColV);
00512     // loads matlab full matrix
00513     static void LoadMatlabTFltVV(const TStr& FNm, TFltVV& MatrixVV);
00514     // prints vector to screen
00515     static void PrintTFltV(const TFltV& Vec, const TStr& VecNm);
00516         // print matrixt to screen
00517         static void PrintTFltVV(const TFltVV& A, const TStr& MatrixNm);
00518     // prints vector to screen
00519     static void PrintTIntV(const TIntV& Vec, const TStr& VecNm);
00520     // fills vector with random numbers
00521     static void FillRnd(TFltV& Vec) { TRnd Rnd(0); FillRnd(Vec, Rnd); }
00522     static void FillRnd(TFltV& Vec, TRnd& Rnd);
00523     // set all components
00524     static void Fill(TFltVV& M, const double& Val);
00525     // sets all compnents to zero
00526     static void FillZero(TFltV& Vec) { Vec.PutAll(0.0); }
00527     static void FillZero(TFltVV& M) { Fill(M, 0.0); }
00528     // set matrix to identity
00529     static void FillIdentity(TFltVV& M);
00530     static void FillIdentity(TFltVV& M, const double& Elt);
00531     // sums elements in vector
00532     static int SumVec(const TIntV& Vec);
00533     static double SumVec(const TFltV& Vec);
00534     // converts full vector to sparse
00535     static void ToSpVec(const TFltV& Vec, TIntFltKdV& SpVec,
00536         const double& CutWordWgtSumPrc = 0.0);
00537     // converts sparse vector to full
00538     static void ToVec(const TIntFltKdV& SpVec, TFltV& Vec, const int& VecLen);
00539 };
00540 
00542 // Template-ised Sparse Operations
00543 template <class TKey, class TDat>
00544 class TSparseOps {
00545 private:
00546         typedef TVec<TKeyDat<TKey, TDat> > TKeyDatV;
00547 public:
00548         static void SparseMerge(const TKeyDatV& SrcV1, const TKeyDatV& SrcV2, TKeyDatV& DstV) {
00549                 DstV.Clr();
00550                 const int Src1Len = SrcV1.Len();
00551                 const int Src2Len = SrcV2.Len();
00552                 int Src1N = 0, Src2N = 0;
00553                 while (Src1N < Src1Len && Src2N < Src2Len) {
00554                         if (SrcV1[Src1N].Key < SrcV2[Src2N].Key) { 
00555                                 DstV.Add(SrcV1[Src1N]); Src1N++;
00556                         } else if (SrcV1[Src1N].Key > SrcV2[Src2N].Key) { 
00557                                 DstV.Add(SrcV2[Src2N]); Src2N++;
00558                         } else { 
00559                                 DstV.Add(TKeyDat<TKey, TDat>(SrcV1[Src1N].Key, SrcV1[Src1N].Dat + SrcV2[Src2N].Dat));
00560                                 Src1N++;  Src2N++; 
00561                         }
00562                 }
00563                 while (Src1N < Src1Len) { DstV.Add(SrcV1[Src1N]); Src1N++; }
00564                 while (Src2N < Src2Len) { DstV.Add(SrcV2[Src2N]); Src2N++; }
00565         }
00566 };
00567 
00568 typedef TSparseOps<TInt, TFlt> TSparseOpsIntFlt;