SNAP Library , Developer Reference  2013-01-07 14:03:36
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     // Result = ||x-y|| (Euclidian)
00199     static double EuclDist(const TFltV& x, const TFltV& y);
00200 
00201     // ||x||^2 (Euclidian)
00202     static double Norm2(const TFltV& x);
00203     // ||x|| (Euclidian)
00204     static double Norm(const TFltV& x);
00205     // x := x / ||x||
00206     static void Normalize(TFltV& x);
00207 
00208     // ||x||^2 (Euclidian), x is sparse
00209     static double Norm2(const TIntFltKdV& x);
00210     // ||x|| (Euclidian), x is sparse
00211     static double Norm(const TIntFltKdV& x);
00212     // x := x / ||x||, x is sparse
00213     static void Normalize(TIntFltKdV& x);
00214 
00215     // ||X(:,ColId)||^2 (Euclidian)
00216     static double Norm2(const TFltVV& X, int ColId);
00217     // ||X(:,ColId)|| (Euclidian)
00218     static double Norm(const TFltVV& X, int ColId);
00219 
00220     // L1 norm of x (Sum[|xi|, i = 1..n])
00221     static double NormL1(const TFltV& x);
00222     // L1 norm of k*x+y (Sum[|k*xi+yi|, i = 1..n])
00223     static double NormL1(double k, const TFltV& x, const TFltV& y);
00224     // L1 norm of x (Sum[|xi|, i = 1..n])
00225     static double NormL1(const TIntFltKdV& x);
00226     // x := x / ||x||_inf
00227     static void NormalizeL1(TFltV& x);
00228     // x := x / ||x||_inf
00229     static void NormalizeL1(TIntFltKdV& x);
00230 
00231     // Linf norm of x (Max{|xi|, i = 1..n})
00232     static double NormLinf(const TFltV& x);
00233     // Linf norm of x (Max{|xi|, i = 1..n})
00234     static double NormLinf(const TIntFltKdV& x);
00235     // x := x / ||x||_inf
00236     static void NormalizeLinf(TFltV& x);
00237     // x := x / ||x||_inf, , x is sparse
00238     static void NormalizeLinf(TIntFltKdV& x);
00239 
00240     // y := k * x
00241     static void MultiplyScalar(const double& k, const TFltV& x, TFltV& y);
00242     // y := k * x
00243     static void MultiplyScalar(const double& k, const TIntFltKdV& x, TIntFltKdV& y);
00244 
00245     // y := A * x
00246     static void Multiply(const TFltVV& A, const TFltV& x, TFltV& y);
00247     // C(:, ColId) := A * x
00248     static void Multiply(const TFltVV& A, const TFltV& x, TFltVV& C, int ColId);
00249     // y := A * B(:, ColId)
00250     static void Multiply(const TFltVV& A, const TFltVV& B, int ColId, TFltV& y);
00251     // C(:, ColIdC) := A * B(:, ColIdB)
00252     static void Multiply(const TFltVV& A, const TFltVV& B, int ColIdB, TFltVV& C, int ColIdC);
00253 
00254     // y := A' * x
00255     static void MultiplyT(const TFltVV& A, const TFltV& x, TFltV& y);
00256 
00257     // C = A * B
00258     static void Multiply(const TFltVV& A, const TFltVV& B, TFltVV& C);
00259         
00260         // D = alpha * A(') * B(') + beta * C(')
00261         typedef enum { GEMM_NO_T = 0, GEMM_A_T = 1, GEMM_B_T = 2, GEMM_C_T = 4 } TLinAlgGemmTranspose;
00262         static void Gemm(const double& Alpha, const TFltVV& A, const TFltVV& B, const double& Beta, 
00263                 const TFltVV& C, TFltVV& D, const int& TransposeFlags);
00264         
00265         // B = A^(-1)
00266         typedef enum { DECOMP_SVD } TLinAlgInverseType;
00267         static void Inverse(const TFltVV& A, TFltVV& B, const TLinAlgInverseType& DecompType);
00268         // subtypes of finding an inverse
00269         static void InverseSVD(const TFltVV& A, TFltVV& B);
00270 
00271         // transpose matrix - B = A'
00272         static void Transpose(const TFltVV& A, TFltVV& B);
00273 
00274     // performes Gram-Schmidt ortogonalization on elements of Q
00275     static void GS(TVec<TFltV>& Q);
00276     // Gram-Schmidt on columns of matrix Q
00277     static void GS(TFltVV& Q);
00278 
00279     // rotates vector (OldX,OldY) for angle Angle (in radians!)
00280     static void Rotate(const double& OldX, const double& OldY, const double& Angle, double& NewX, double& NewY);
00281 
00282     // checks if set of vectors is ortogonal
00283     static void AssertOrtogonality(const TVec<TFltV>& Vecs, const double& Threshold);
00284     static void AssertOrtogonality(const TFltVV& Vecs, const double& Threshold);
00285 };
00286 
00288 // Numerical-Recepies-Exception
00289 class TNSException {
00290 public:
00291     TStr Message;
00292 public:
00293     TNSException(const TStr& Msg): Message(Msg) {}
00294 };
00295 
00297 // Numerical-Linear-Algebra (copied from Numerical Recepies)
00298 class TNumericalStuff {
00299 private:
00300   static double sqr(double a);
00301   static double sign(double a, double b);
00302 
00303   // Computes (a^2 + b^2)^(1/2) without
00304   // destructive underflow or overflow.
00305   static double pythag(double a, double b);
00306 
00307   //displays error message to screen
00308   static void nrerror(const TStr& error_text);
00309 
00310 public:
00311     // Householder reduction of a real, symmetric matrix a[1..n][1..n].
00312     // On output, a is replaced by the orthogonal matrix Q eecting the
00313     // transformation. d[1..n] returns the diagonal elements of the
00314     // tridiagonal matrix, and e[1..n] the o-diagonal elements, with
00315     // e[1]=0. Several statements, as noted in comments, can be omitted
00316     // if only eigenvalues are to be found, in which case a contains no
00317     // useful information on output. Otherwise they are to be included.
00318     static void SymetricToTridiag(TFltVV& a, int n, TFltV& d, TFltV& e);
00319 
00320         // QL algorithm with implicit shifts, to determine the eigenvalues
00321         // and eigenvectors of a real, symmetric, tridiagonal matrix, or of
00322         // a real, symmetric matrix previously reduced by tred2 x11.2. On
00323         // input, d[1..n] contains the diagonal elements of the tridiagonal
00324         // matrix. On output, it returns the eigenvalues. The vector e[1..n]
00325         // inputs the subdiagonal elements of the tridiagonal matrix, with
00326         // e[1] arbitrary. On output e is destroyed. When finding only the
00327         // eigenvalues, several lines may be omitted, as noted in the comments.
00328         // If the eigenvectors of a tridiagonal matrix are desired, the matrix
00329         // z[1..n][1..n] is input as the identity matrix. If the eigenvectors
00330         // of a matrix that has been reduced by tred2 are required, then z is
00331         // input as the matrix output by tred2. In either case, the kth column
00332         // of z returns the normalized eigenvector corresponding to d[k].
00333         static void EigSymmetricTridiag(TFltV& d, TFltV& e, int n, TFltVV& z);
00334 
00335         // Given a positive-dedinite symmetric matrix A(n,n), this routine
00336         // constructs its Cholesky decomposition, A = L * L^T . On input, only
00337         // the upper triangle of A need be given; it is not modified. The
00338         // Cholesky factor L is returned in the lower triangle of A, except for
00339         // its diagonal elements which are returned in p(n).
00340         static void CholeskyDecomposition(TFltVV& A, TFltV& p);
00341 
00342         // Solves the set of n linear equations A * x = b, where A is a
00343         // positive-definite symmetric matrix. A(n,n) and p[1..n] are input
00344         // as the output of the routine choldc. Only the lower triangle of A
00345         // is accessed. b(n) is input as the right-hand side vector. The
00346         // solution vector is returned in x(n). A  and p are not modified and
00347         // can be left in place for successive calls with diferent right-hand
00348         // sides b. b is not modified unless you identify b and x in the calling
00349         // sequence, which is allowed.
00350         static void CholeskySolve(const TFltVV& A, const TFltV& p, const TFltV& b, TFltV& x);
00351 
00352         // Solves system of linear equations A * x = b, where A is symetric
00353         // positive-definite matrix. A is first decomposed using
00354         // CholeskyDecomposition and after solved using CholeskySolve. Only
00355         // upper triangle of A need be given and it is not modified. However,
00356         // lower triangle is modified!
00357         static void SolveSymetricSystem(TFltVV& A, const TFltV& b, TFltV& x);
00358 
00359     // solve system A x_i = e_i for i = 1..n, where A and p are output
00360     // from CholeskyDecomposition. Result is stored to upper triangule
00361     // (possible since inverse of symetric matrix is also symetric! Sigh...)
00362     static void InverseSubstitute(TFltVV& A, const TFltV& p);
00363 
00364     // Calculates inverse of symetric positiv definit matrix
00365     // Matrix is given as upper triangule of A, result is stored
00366     // in upper triangule of A. Lower triangule is random (actually
00367     // it has part of Choleksy decompositon of A)
00368     static void InverseSymetric(TFltVV& A);
00369 
00370     // calcualtes inverse of upper triagonal matrix A
00371     // lower triangle is messed up...
00372     static void InverseTriagonal(TFltVV& A);
00373 
00374     // Given a matrix a[1..n][1..n], this routine replaces it by the LU
00375     // decomposition of a rowwise permutation of itself. a and n are input.
00376     // a is output, arranged as in equation (2.3.14) above; indx[1..n] is
00377     // an output vector that records the row permutation efected by the partial
00378     // pivoting; d is output as +-1 depending on whether the number of row
00379     // interchanges was even or odd, respectively. This routine is used in
00380     // combination with lubksb to solve linear equations or invert a matrix.
00381     static void LUDecomposition(TFltVV& A, TIntV& indx, double& d);
00382 
00383     // Solves the set of n linear equations A*X = B. Here a[1..n][1..n] is input,
00384     // not as the matrix A but rather as its LU decomposition, determined by the
00385     // routine ludcmp. indx[1..n] is input as the permutation vector returned by
00386     // ludcmp. b[1..n] is input as the right-hand side vector B, and returns with
00387     // the solution vector X. a, n, and indx are not modified by this routine and
00388     // can be left in place for successive calls with diferent right-hand sides b.
00389     // This routine takes into account the possibility that b will begin with many
00390     // zero elements, so it is efficient for use in matrix inversion.
00391     static void LUSolve(const TFltVV& A, const TIntV& indx, TFltV& b);
00392 
00393     // Solves system of linear equations A * x = b. A is first decomposed using
00394     // LUDecomposition and after solved using LUSolve. A is modified!
00395     static void SolveLinearSystem(TFltVV& A, const TFltV& b, TFltV& x);
00396 };
00397 
00399 // Sparse-SVD
00400 //   Calculates singular-value-decompositon for sparse matrixes.
00401 //   If A is a matrix than A is decomposed to A = U S V'
00402 //   where S is diagonal with singular values on diagonal and U
00403 //   and V are ortogonal (U'*U = V'*V = I).
00404 typedef enum { ssotNoOrto, ssotSelective, ssotFull } TSpSVDReOrtoType;
00405 class TSparseSVD {
00406 private:
00407     // Result = Matrix' * Matrix * Vec(:,ColId)
00408     static void MultiplyATA(const TMatrix& Matrix,
00409         const TFltVV& Vec, int ColId, TFltV& Result);
00410     // Result = Matrix' * Matrix * Vec
00411     static void MultiplyATA(const TMatrix& Matrix,
00412         const TFltV& Vec, TFltV& Result);
00413 public:
00414     // calculates NumEig eigen values of symetric matrix
00415     // if SvdMatrixProductP than matrix Matrix'*Matrix is used
00416     static void SimpleLanczos(const TMatrix& Matrix,
00417         const int& NumEig, TFltV& EigValV,
00418         const bool& DoLocalReortoP = false,
00419         const bool& SvdMatrixProductP = false);
00420     // fast, calculates NumEig largers eigen values and vectors
00421     // kk should be something like 4*NumEig
00422     // if SvdMatrixProductP than matrix Matrix'*Matrix is used
00423     static void Lanczos(const TMatrix& Matrix,
00424         int NumEig, int Iters, const TSpSVDReOrtoType& ReOrtoType,
00425         TFltV& EigValV, TFltVV& EigVecVV,
00426         const bool& SvdMatrixProductP = false);
00427     static void Lanczos2(const TMatrix& Matrix,
00428         int MaxNumEig, int MaxSecs, const TSpSVDReOrtoType& ReOrtoType,
00429         TFltV& EigValV, TFltVV& EigVecVV,
00430         const bool& SvdMatrixProductP = false);
00431 
00432     // calculates only singular values (based on SimpleLanczos)
00433     static void SimpleLanczosSVD(const TMatrix& Matrix,
00434         const int& CalcSV, TFltV& SngValV,
00435         const bool& DoLocalReortoP = false);
00436     // fast, calculates NumSV largers SV (based on Lanczos)
00437     static void LanczosSVD(const TMatrix& Matrix,
00438         int NumSV, int Iters, const TSpSVDReOrtoType& ReOrtoType,
00439         TFltV& SgnValV, TFltVV& LeftSgnVecVV, TFltVV& RightSgnVecVV);
00440 
00441     // slow - ortogonal iteration
00442     static void OrtoIterSVD(const TMatrix& Matrix,
00443         int NumSV, int IterN, TFltV& SgnValV);
00444 
00445     // projects sparse vector to space spanned by columns of matrix U
00446     static void Project(const TIntFltKdV& Vec, const TFltVV& U, TFltV& ProjVec);
00447 };
00448 
00450 // Sigmoid  --  made by Janez(TM)
00451 //  (y = 1/[1 + exp[-Ax+B]])
00452 class TSigmoid {
00453 private:
00454     TFlt A;
00455     TFlt B;
00456 private:
00457   // Evaluates how well the sigmoid function fits the data.
00458   // J(A, B) = - ln prod_i P(Y = y_i | Z = z_i).  The 'data' parameter
00459   // should contain (z_i, y_i) pairs.  Smaller J means a better fit.
00460   static double EvaluateFit(const TFltIntKdV& data, const double A, const double B);
00461   // Computes not only J but also its partial derivatives.
00462   static void EvaluateFit(const TFltIntKdV& data, const double A,
00463         const double B, double& J, double& JA, double& JB);
00464   // Let J(lambda) = J(A + lambda U, B + lambda V).
00465     // This function computes J and its first and second derivatives.
00466   // They can be used to choose a good lambda (using Newton's method)
00467     // when minimizing J. -- This method has not been tested yet.
00468   static void EvaluateFit(const TFltIntKdV& data, const double A,
00469         const double B, const double U, const double V, const double lambda,
00470     double& J, double& JJ, double& JJJ);
00471 public:
00472     TSigmoid() { };
00473     TSigmoid(const double& A_, const double& B_): A(A_), B(B_) { };
00474         // Tries to find a pair (A, B) that minimizes J(A, B).
00475     // Uses gradient descent.
00476     TSigmoid(const TFltIntKdV& data);
00477 
00478     TSigmoid(TSIn& SIn) { A.Load(SIn); B.Load(SIn); }
00479     void Load(TSIn& SIn) { A.Load(SIn); B.Load(SIn); }
00480     void Save(TSOut& SOut) const {A.Save(SOut); B.Save(SOut);}
00481 
00482     double GetVal(const double& x) const {
00483         return 1.0 / (1.0 + exp(-A * x + B)); }
00484     double operator()(const double& x) const {
00485         return GetVal(x); }
00486 
00487     void GetSigmoidAB(double& A_, double& B_) { A_=A; B_=B; }
00488 };
00489 
00491 // Useful stuff (hopefuly)
00492 class TLAMisc {
00493 public:
00494         // Dumps vector to file so Excel can read it
00495     static void SaveCsvTFltV(const TFltV& Vec, TSOut& SOut);
00496         // Dumps sparse vector to file so Matlab can read it
00497     static void SaveMatlabTFltIntKdV(const TIntFltKdV& SpV, const int& ColN, TSOut& SOut);
00498         // Dumps vector to file so Matlab can read it
00499     static void SaveMatlabTFltV(const TFltV& m, const TStr& FName);
00500         // Dumps vector to file so Matlab can read it
00501     static void SaveMatlabTIntV(const TIntV& m, const TStr& FName);
00502         // Dumps column ColId from m to file so Matlab can read it
00503     static void SaveMatlabTFltVVCol(const TFltVV& m, int ColId, const TStr& FName);
00504         // Dumps matrix to file so Matlab can read it
00505     static void SaveMatlabTFltVV(const TFltVV& m, const TStr& FName);
00506         // Dumps main minor rowN x colN to file so Matlab can read it
00507         static void SaveMatlabTFltVVMjrSubMtrx(const TFltVV& m, int rowN, int colN, const TStr& FName);
00508     // loads matlab full matrix
00509     static void LoadMatlabTFltVV(const TStr& FNm, TVec<TFltV>& ColV);
00510     // loads matlab full matrix
00511     static void LoadMatlabTFltVV(const TStr& FNm, TFltVV& MatrixVV);
00512     // prints vector to screen
00513     static void PrintTFltV(const TFltV& Vec, const TStr& VecNm);
00514         // print matrixt to screen
00515         static void PrintTFltVV(const TFltVV& A, const TStr& MatrixNm);
00516     // prints vector to screen
00517     static void PrintTIntV(const TIntV& Vec, const TStr& VecNm);
00518     // fills vector with random numbers
00519     static void FillRnd(TFltV& Vec) { TRnd Rnd(0); FillRnd(Vec, Rnd); }
00520     static void FillRnd(TFltV& Vec, TRnd& Rnd);
00521     // set all components
00522     static void Fill(TFltVV& M, const double& Val);
00523     // sets all compnents to zero
00524     static void FillZero(TFltV& Vec) { Vec.PutAll(0.0); }
00525     static void FillZero(TFltVV& M) { Fill(M, 0.0); }
00526     // set matrix to identity
00527     static void FillIdentity(TFltVV& M);
00528     static void FillIdentity(TFltVV& M, const double& Elt);
00529     // sums elements in vector
00530     static int SumVec(const TIntV& Vec);
00531     static double SumVec(const TFltV& Vec);
00532     // converts full vector to sparse
00533     static void ToSpVec(const TFltV& Vec, TIntFltKdV& SpVec,
00534         const double& CutWordWgtSumPrc = 0.0);
00535     // converts sparse vector to full
00536     static void ToVec(const TIntFltKdV& SpVec, TFltV& Vec, const int& VecLen);
00537 };
00538 
00540 // Template-ised Sparse Operations
00541 template <class TKey, class TDat>
00542 class TSparseOps {
00543 private:
00544         typedef TVec<TKeyDat<TKey, TDat> > TKeyDatV;
00545 public:
00546         static void SparseMerge(const TKeyDatV& SrcV1, const TKeyDatV& SrcV2, TKeyDatV& DstV) {
00547                 DstV.Clr();
00548                 const int Src1Len = SrcV1.Len();
00549                 const int Src2Len = SrcV2.Len();
00550                 int Src1N = 0, Src2N = 0;
00551                 while (Src1N < Src1Len && Src2N < Src2Len) {
00552                         if (SrcV1[Src1N].Key < SrcV2[Src2N].Key) { 
00553                                 DstV.Add(SrcV1[Src1N]); Src1N++;
00554                         } else if (SrcV1[Src1N].Key > SrcV2[Src2N].Key) { 
00555                                 DstV.Add(SrcV2[Src2N]); Src2N++;
00556                         } else { 
00557                                 DstV.Add(TKeyDat<TKey, TDat>(SrcV1[Src1N].Key, SrcV1[Src1N].Dat + SrcV2[Src2N].Dat));
00558                                 Src1N++;  Src2N++; 
00559                         }
00560                 }
00561                 while (Src1N < Src1Len) { DstV.Add(SrcV1[Src1N]); Src1N++; }
00562                 while (Src2N < Src2Len) { DstV.Add(SrcV2[Src2N]); Src2N++; }
00563         }
00564 };
00565 
00566 typedef TSparseOps<TInt, TFlt> TSparseOpsIntFlt;