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
TLinAlg Class Reference

#include <linalg.h>

List of all members.

Public Types

enum  TLinAlgGemmTranspose { GEMM_NO_T = 0, GEMM_A_T = 1, GEMM_B_T = 2, GEMM_C_T = 4 }
enum  TLinAlgInverseType { DECOMP_SVD }

Static Public Member Functions

static double DotProduct (const TFltV &x, const TFltV &y)
static double DotProduct (const TFltVV &X, int ColIdX, const TFltVV &Y, int ColIdY)
static double DotProduct (const TFltVV &X, int ColId, const TFltV &Vec)
static double DotProduct (const TIntFltKdV &x, const TIntFltKdV &y)
static double DotProduct (const TFltV &x, const TIntFltKdV &y)
static double DotProduct (const TFltVV &X, int ColId, const TIntFltKdV &y)
static void LinComb (const double &p, const TFltV &x, const double &q, const TFltV &y, TFltV &z)
static void ConvexComb (const double &p, const TFltV &x, const TFltV &y, TFltV &z)
static void AddVec (const double &k, const TFltV &x, const TFltV &y, TFltV &z)
static void AddVec (const double &k, const TIntFltKdV &x, const TFltV &y, TFltV &z)
static void AddVec (const double &k, const TIntFltKdV &x, TFltV &y)
static void AddVec (double k, const TFltVV &X, int ColIdX, TFltVV &Y, int ColIdY)
static void AddVec (double k, const TFltVV &X, int ColId, TFltV &Result)
static void AddVec (const TIntFltKdV &x, const TIntFltKdV &y, TIntFltKdV &z)
static double SumVec (const TFltV &x)
static double SumVec (double k, const TFltV &x, const TFltV &y)
static double EuclDist2 (const TFltV &x, const TFltV &y)
static double EuclDist2 (const TFltPr &x, const TFltPr &y)
static double EuclDist (const TFltV &x, const TFltV &y)
static double EuclDist (const TFltPr &x, const TFltPr &y)
static double Norm2 (const TFltV &x)
static double Norm (const TFltV &x)
static void Normalize (TFltV &x)
static double Norm2 (const TIntFltKdV &x)
static double Norm (const TIntFltKdV &x)
static void Normalize (TIntFltKdV &x)
static double Norm2 (const TFltVV &X, int ColId)
static double Norm (const TFltVV &X, int ColId)
static double NormL1 (const TFltV &x)
static double NormL1 (double k, const TFltV &x, const TFltV &y)
static double NormL1 (const TIntFltKdV &x)
static void NormalizeL1 (TFltV &x)
static void NormalizeL1 (TIntFltKdV &x)
static double NormLinf (const TFltV &x)
static double NormLinf (const TIntFltKdV &x)
static void NormalizeLinf (TFltV &x)
static void NormalizeLinf (TIntFltKdV &x)
static void MultiplyScalar (const double &k, const TFltV &x, TFltV &y)
static void MultiplyScalar (const double &k, const TIntFltKdV &x, TIntFltKdV &y)
static void Multiply (const TFltVV &A, const TFltV &x, TFltV &y)
static void Multiply (const TFltVV &A, const TFltV &x, TFltVV &C, int ColId)
static void Multiply (const TFltVV &A, const TFltVV &B, int ColId, TFltV &y)
static void Multiply (const TFltVV &A, const TFltVV &B, int ColIdB, TFltVV &C, int ColIdC)
static void MultiplyT (const TFltVV &A, const TFltV &x, TFltV &y)
static void Multiply (const TFltVV &A, const TFltVV &B, TFltVV &C)
static void Gemm (const double &Alpha, const TFltVV &A, const TFltVV &B, const double &Beta, const TFltVV &C, TFltVV &D, const int &TransposeFlags)
static void Inverse (const TFltVV &A, TFltVV &B, const TLinAlgInverseType &DecompType)
static void InverseSVD (const TFltVV &A, TFltVV &B)
static void Transpose (const TFltVV &A, TFltVV &B)
static void GS (TVec< TFltV > &Q)
static void GS (TFltVV &Q)
static void Rotate (const double &OldX, const double &OldY, const double &Angle, double &NewX, double &NewY)
static void AssertOrtogonality (const TVec< TFltV > &Vecs, const double &Threshold)
static void AssertOrtogonality (const TFltVV &Vecs, const double &Threshold)

Detailed Description

Definition at line 156 of file linalg.h.


Member Enumeration Documentation

Enumerator:
GEMM_NO_T 
GEMM_A_T 
GEMM_B_T 
GEMM_C_T 

Definition at line 263 of file linalg.h.

Enumerator:
DECOMP_SVD 

Definition at line 268 of file linalg.h.


Member Function Documentation

void TLinAlg::AddVec ( const double &  k,
const TFltV x,
const TFltV y,
TFltV z 
) [static]

Definition at line 232 of file linalg.cpp.

References LinComb().

Referenced by GS(), TSparseSVD::Lanczos(), TSparseSVD::Lanczos2(), TFullColMatrix::PMultiply(), and TSparseSVD::SimpleLanczos().

                                                                              {
    LinComb(k, x, 1.0, y, z);
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::AddVec ( const double &  k,
const TIntFltKdV x,
const TFltV y,
TFltV z 
) [static]

Definition at line 236 of file linalg.cpp.

References Assert, and TVec< TVal, TSizeTy >::Len().

                                                                                   {
    Assert(y.Len() == z.Len());
    z = y; // first we set z to be y
    // and than we add x to z (==y)
    const int xLen = x.Len(), yLen = y.Len();
    for (int i = 0; i < xLen; i++) {
        const int ii = x[i].Key;
        if (ii < yLen) {
            z[ii] = k * x[i].Dat + y[ii];
        }
    }
}

Here is the call graph for this function:

void TLinAlg::AddVec ( const double &  k,
const TIntFltKdV x,
TFltV y 
) [static]

Definition at line 249 of file linalg.cpp.

References TVec< TVal, TSizeTy >::Len().

                                                                   {
    const int xLen = x.Len(), yLen = y.Len();
    for (int i = 0; i < xLen; i++) {
        const int ii = x[i].Key;
        if (ii < yLen) {
            y[ii] += k * x[i].Dat;
        }
    }
}

Here is the call graph for this function:

void TLinAlg::AddVec ( double  k,
const TFltVV X,
int  ColIdX,
TFltVV Y,
int  ColIdY 
) [static]

Definition at line 259 of file linalg.cpp.

References Assert, and TVVec< TVal >::GetRows().

                                                                                {
    Assert(X.GetRows() == Y.GetRows());
    const int len = Y.GetRows();
    for (int i = 0; i < len; i++) {
        Y(i,ColIdY) = Y(i,ColIdY) + k * X(i, ColIdX);
    }
}

Here is the call graph for this function:

void TLinAlg::AddVec ( double  k,
const TFltVV X,
int  ColId,
TFltV Result 
) [static]

Definition at line 267 of file linalg.cpp.

References Assert, TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

                                                                       {
    Assert(X.GetRows() == Result.Len());
    const int len = Result.Len();
    for (int i = 0; i < len; i++) {
        Result[i] = Result[i] + k * X(i, ColId);
    }
}

Here is the call graph for this function:

void TLinAlg::AddVec ( const TIntFltKdV x,
const TIntFltKdV y,
TIntFltKdV z 
) [static]

Definition at line 275 of file linalg.cpp.

References TSparseOps< TKey, TDat >::SparseMerge().

Here is the call graph for this function:

void TLinAlg::AssertOrtogonality ( const TVec< TFltV > &  Vecs,
const double &  Threshold 
) [static]

Definition at line 616 of file linalg.cpp.

References TFlt::Abs(), DotProduct(), TVec< TVal, TSizeTy >::Len(), and Norm2().

                                                                                 {
    int m = Vecs.Len();
    for (int i = 0; i < m; i++) {
        for (int j = 0; j < i; j++) {
            double res = DotProduct(Vecs[i], Vecs[j]);
            if (TFlt::Abs(res) > Threshold)
                printf("<%d,%d> = %.5f", i,j,res);
        }
        double norm = Norm2(Vecs[i]);
        if (TFlt::Abs(norm-1) > Threshold)
            printf("||%d|| = %.5f", i, norm);
    }
}

Here is the call graph for this function:

void TLinAlg::AssertOrtogonality ( const TFltVV Vecs,
const double &  Threshold 
) [static]

Definition at line 630 of file linalg.cpp.

References TFlt::Abs(), DotProduct(), TVVec< TVal >::GetCols(), and Norm2().

                                                                            {
    int m = Vecs.GetCols();
    for (int i = 0; i < m; i++) {
        for (int j = 0; j < i; j++) {
            double res = DotProduct(Vecs, i, Vecs, j);
            if (TFlt::Abs(res) > Threshold)
                printf("<%d,%d> = %.5f", i, j, res);
        }
        double norm = Norm2(Vecs, i);
        if (TFlt::Abs(norm-1) > Threshold)
            printf("||%d|| = %.5f", i, norm);
    }
    printf("\n");
}

Here is the call graph for this function:

void TLinAlg::ConvexComb ( const double &  p,
const TFltV x,
const TFltV y,
TFltV z 
) [static]

Definition at line 227 of file linalg.cpp.

References AssertR, TFlt::GetStr(), and LinComb().

                                                                                  {
    AssertR(0.0 <= p && p <= 1.0, TFlt::GetStr(p));
    LinComb(p, x, 1.0 - p, y, z);
}

Here is the call graph for this function:

double TLinAlg::DotProduct ( const TFltV x,
const TFltV y 
) [static]

Definition at line 165 of file linalg.cpp.

References TStr::Fmt(), IAssertR, and TVec< TVal, TSizeTy >::Len().

Referenced by AssertOrtogonality(), TLogRegFit::GetStepSizeByLineSearch(), TAGMFit::GetStepSizeByLineSearchForLambda(), GS(), TSparseSVD::Lanczos(), TSparseSVD::Lanczos2(), TLogRegFit::MLENewton(), Norm2(), TFullColMatrix::PMultiplyT(), and TSparseSVD::SimpleLanczos().

                                                         {
    IAssertR(x.Len() == y.Len(), TStr::Fmt("%d != %d", x.Len(), y.Len()));
    double result = 0.0; int Len = x.Len();
    for (int i = 0; i < Len; i++)
        result += x[i] * y[i];
    return result;
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::DotProduct ( const TFltVV X,
int  ColIdX,
const TFltVV Y,
int  ColIdY 
) [static]

Definition at line 173 of file linalg.cpp.

References Assert, and TVVec< TVal >::GetRows().

                                                                                   {
    Assert(X.GetRows() == Y.GetRows());
    double result = 0.0, len = X.GetRows();
    for (int i = 0; i < len; i++)
        result = result + X(i,ColIdX) * Y(i,ColIdY);
    return result;
}

Here is the call graph for this function:

double TLinAlg::DotProduct ( const TFltVV X,
int  ColId,
const TFltV Vec 
) [static]

Definition at line 181 of file linalg.cpp.

References Assert, TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

                                                                       {
    Assert(X.GetRows() == Vec.Len());
    double result = 0.0; int Len = X.GetRows();
    for (int i = 0; i < Len; i++)
        result += X(i,ColId) * Vec[i];
    return result;
}

Here is the call graph for this function:

double TLinAlg::DotProduct ( const TIntFltKdV x,
const TIntFltKdV y 
) [static]

Definition at line 189 of file linalg.cpp.

References TVec< TVal, TSizeTy >::Len().

                                                                   {
    const int xLen = x.Len(), yLen = y.Len();
    double Res = 0.0; int i1 = 0, i2 = 0;
    while (i1 < xLen && i2 < yLen) {
        if (x[i1].Key < y[i2].Key) i1++;
        else if (x[i1].Key > y[i2].Key) i2++;
        else { Res += x[i1].Dat * y[i2].Dat;  i1++;  i2++; }
    }
    return Res;
}

Here is the call graph for this function:

double TLinAlg::DotProduct ( const TFltV x,
const TIntFltKdV y 
) [static]

Definition at line 200 of file linalg.cpp.

References TVec< TVal, TSizeTy >::Len().

                                                              {
    double Res = 0.0; const int xLen = x.Len(), yLen = y.Len();
    for (int i = 0; i < yLen; i++) {
        const int key = y[i].Key;
        if (key < xLen) Res += y[i].Dat * x[key];
    }
    return Res;
}

Here is the call graph for this function:

double TLinAlg::DotProduct ( const TFltVV X,
int  ColId,
const TIntFltKdV y 
) [static]

Definition at line 209 of file linalg.cpp.

References TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

                                                                          {
    double Res = 0.0; const int n = X.GetRows(), yLen = y.Len();
    for (int i = 0; i < yLen; i++) {
        const int key = y[i].Key;
        if (key < n) Res += y[i].Dat * X(key,ColId);
    }
    return Res;
}

Here is the call graph for this function:

double TLinAlg::EuclDist ( const TFltV x,
const TFltV y 
) [static]

Definition at line 312 of file linalg.cpp.

References EuclDist2().

                                                       {
    return sqrt(EuclDist2(x, y));
}

Here is the call graph for this function:

double TLinAlg::EuclDist ( const TFltPr x,
const TFltPr y 
) [static]

Definition at line 316 of file linalg.cpp.

References EuclDist2().

                                                         {
    return sqrt(EuclDist2(x, y));
}

Here is the call graph for this function:

double TLinAlg::EuclDist2 ( const TFltV x,
const TFltV y 
) [static]

Definition at line 298 of file linalg.cpp.

References Assert, TVec< TVal, TSizeTy >::Len(), and TMath::Sqr().

Referenced by EuclDist().

                                                        {
    Assert(x.Len() == y.Len());
    const int len = x.Len();
    double Res = 0.0;
    for (int i = 0; i < len; i++) {
        Res += TMath::Sqr(x[i] - y[i]);
    }
    return Res;
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::EuclDist2 ( const TFltPr x,
const TFltPr y 
) [static]

Definition at line 308 of file linalg.cpp.

References TMath::Sqr(), TPair< TVal1, TVal2 >::Val1, and TPair< TVal1, TVal2 >::Val2.

                                                          {
    return TMath::Sqr(x.Val1 - y.Val1) + TMath::Sqr(x.Val2 - y.Val2);
}

Here is the call graph for this function:

void TLinAlg::Gemm ( const double &  Alpha,
const TFltVV A,
const TFltVV B,
const double &  Beta,
const TFltVV C,
TFltVV D,
const int &  TransposeFlags 
) [static]

Definition at line 493 of file linalg.cpp.

References Assert, TVVec< TVal >::At(), GEMM_A_T, GEMM_B_T, GEMM_C_T, TVVec< TVal >::GetCols(), and TVVec< TVal >::GetRows().

                                                                       {

        bool tA = (TransposeFlags & GEMM_A_T) == GEMM_A_T;
        bool tB = (TransposeFlags & GEMM_B_T) == GEMM_B_T;
        bool tC = (TransposeFlags & GEMM_C_T) == GEMM_C_T;

        // setting dimensions
        int a_i = tA ? A.GetRows() : A.GetCols();
        int a_j = tA ? A.GetCols() : A.GetRows();

        int b_i = tB ? A.GetRows() : A.GetCols();
        int b_j = tB ? A.GetCols() : A.GetRows();

        int c_i = tC ? A.GetRows() : A.GetCols();
        int c_j = tC ? A.GetCols() : A.GetRows();

        int d_i = D.GetCols();
        int d_j = D.GetRows();
        
        // assertions for dimensions
  bool Cnd = a_i == c_j && b_i == c_i && a_i == b_j && c_i == d_i && c_j == d_j;
  if (Cnd) {
          Assert(Cnd);
  }

        double Aij, Bij, Cij;

        // rows
        for (int j = 0; j < a_j; j++) {
                // cols
                for (int i = 0; i < a_i; i++) {
                        // not optimized for speed (!)
                        double sum = 0;
                        for (int k = 0; k < a_i; k++) {
                                Aij = tA ? A.At(j, k) : A.At(k, j);
                                Bij = tB ? B.At(k, i) : B.At(i, k);
                                sum += Alpha * Aij * Bij;
                        }
                        Cij = tC ? C.At(i, j) : C.At(j, i);
                        sum += Beta * Cij;
                        D.At(i, j) = sum;
                }
        }
}

Here is the call graph for this function:

void TLinAlg::GS ( TVec< TFltV > &  Q) [static]

Definition at line 584 of file linalg.cpp.

References AddVec(), DotProduct(), IAssert, TVec< TVal, TSizeTy >::Len(), and Normalize().

Referenced by TSparseSVD::OrtoIterSVD().

                               {
    IAssert(Q.Len() > 0);
    int m = Q.Len(); // int n = Q[0].Len();
    for (int i = 0; i < m; i++) {
        printf("%d\r",i);
        for (int j = 0; j < i; j++) {
            double r = TLinAlg::DotProduct(Q[i], Q[j]);
            TLinAlg::AddVec(-r,Q[j],Q[i],Q[i]);
        }
        TLinAlg::Normalize(Q[i]);
    }
    printf("\n");
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::GS ( TFltVV Q) [static]

Definition at line 598 of file linalg.cpp.

References AddVec(), DotProduct(), TVVec< TVal >::GetCols(), TVVec< TVal >::GetRows(), and Norm().

                          {
    int m = Q.GetCols(), n = Q.GetRows();
    for (int i = 0; i < m; i++) {
        for (int j = 0; j < i; j++) {
            double r = TLinAlg::DotProduct(Q, i, Q, j);
            TLinAlg::AddVec(-r,Q,j,Q,i);
        }
        double nr = TLinAlg::Norm(Q,i);
        for (int k = 0; k < n; k++)
            Q(k,i) = Q(k,i) / nr;
    }
}

Here is the call graph for this function:

void TLinAlg::Inverse ( const TFltVV A,
TFltVV B,
const TLinAlgInverseType DecompType 
) [static]

Definition at line 548 of file linalg.cpp.

References DECOMP_SVD, and InverseSVD().

                                                                                      {
        switch (DecompType) {
                case DECOMP_SVD:
                        InverseSVD(A, B);
        }
}

Here is the call graph for this function:

void TLinAlg::InverseSVD ( const TFltVV A,
TFltVV B 
) [static]

Definition at line 555 of file linalg.cpp.

References TVVec< TVal >::At(), TVVec< TVal >::Gen(), TVVec< TVal >::GetCols(), TVVec< TVal >::GetRows(), TVec< TVal, TSizeTy >::Len(), and TSvd::Svd().

Referenced by Inverse().

                                                   {
        // create temp matrices
        TFltVV U, V;
        TFltV E;
        TSvd SVD;

        U.Gen(M.GetRows(), M.GetRows());        
        V.Gen(M.GetCols(), M.GetCols());

        // do the SVD decompostion
        SVD.Svd(M, U, E, V);

        // calculate reciprocal values for diagonal matrix = inverse diagonal
        for (int i = 0; i < E.Len(); i++) {
                E[i] = 1 / E[i];
        }

        // calculate pseudoinverse: M^(-1) = V * E^(-1) * U'
        for (int i = 0; i < U.GetCols(); i++) {
                for (int j = 0; j < V.GetRows(); j++) {
                        double sum = 0;
                        for (int k = 0; k < U.GetCols(); k++) {
                                sum += E[k] * V.At(i, k) * U.At(j, k);
                        }
                        B.At(i, j) = sum;
                }
        }       
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::LinComb ( const double &  p,
const TFltV x,
const double &  q,
const TFltV y,
TFltV z 
) [static]

Definition at line 218 of file linalg.cpp.

References Assert, and TVec< TVal, TSizeTy >::Len().

Referenced by AddVec(), and ConvexComb().

                                                   {

    Assert(x.Len() == y.Len() && y.Len() == z.Len());
    const int Len = x.Len();
    for (int i = 0; i < Len; i++) {
        z[i] = p * x[i] + q * y[i]; }
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::Multiply ( const TFltVV A,
const TFltV x,
TFltV y 
) [static]

Definition at line 428 of file linalg.cpp.

References Assert, TVVec< TVal >::GetCols(), TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

Referenced by TSparseSVD::Lanczos(), and TSparseSVD::Lanczos2().

                                                                {
    Assert(A.GetCols() == x.Len() && A.GetRows() == y.Len());
    int n = A.GetRows(), m = A.GetCols();
    for (int i = 0; i < n; i++) {
        y[i] = 0.0;
        for (int j = 0; j < m; j++)
            y[i] += A(i,j) * x[j];
    }
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::Multiply ( const TFltVV A,
const TFltV x,
TFltVV C,
int  ColId 
) [static]

Definition at line 438 of file linalg.cpp.

References Assert, TVVec< TVal >::GetCols(), TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

                                                                            {
    Assert(A.GetCols() == x.Len() && A.GetRows() == C.GetRows());
    int n = A.GetRows(), m = A.GetCols();
    for (int i = 0; i < n; i++) {
        C(i,ColId) = 0.0;
        for (int j = 0; j < m; j++)
            C(i,ColId) += A(i,j) * x[j];
    }
}

Here is the call graph for this function:

void TLinAlg::Multiply ( const TFltVV A,
const TFltVV B,
int  ColId,
TFltV y 
) [static]

Definition at line 448 of file linalg.cpp.

References Assert, TVVec< TVal >::GetCols(), TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

                                                                            {
    Assert(A.GetCols() == B.GetRows() && A.GetRows() == y.Len());
    int n = A.GetRows(), m = A.GetCols();
    for (int i = 0; i < n; i++) {
        y[i] = 0.0;
        for (int j = 0; j < m; j++)
            y[i] += A(i,j) * B(j,ColId);
    }
}

Here is the call graph for this function:

void TLinAlg::Multiply ( const TFltVV A,
const TFltVV B,
int  ColIdB,
TFltVV C,
int  ColIdC 
) [static]

Definition at line 458 of file linalg.cpp.

References Assert, TVVec< TVal >::GetCols(), and TVVec< TVal >::GetRows().

                                                                                         {
    Assert(A.GetCols() == B.GetRows() && A.GetRows() == C.GetRows());
    int n = A.GetRows(), m = A.GetCols();
    for (int i = 0; i < n; i++) {
        C(i,ColIdC) = 0.0;
        for (int j = 0; j < m; j++)
            C(i,ColIdC) += A(i,j) * B(j,ColIdB);
    }
}

Here is the call graph for this function:

void TLinAlg::Multiply ( const TFltVV A,
const TFltVV B,
TFltVV C 
) [static]

Definition at line 478 of file linalg.cpp.

References Assert, TVVec< TVal >::GetCols(), and TVVec< TVal >::GetRows().

                                                                  {
    Assert(A.GetRows() == C.GetRows() && B.GetCols() == C.GetCols() && A.GetCols() == B.GetRows());

    int n = C.GetRows(), m = C.GetCols(), l = A.GetCols();
    for (int i = 0; i < n; i++) {
        for (int j = 0; j < m; j++) {
            double sum = 0.0;
            for (int k = 0; k < l; k++)
                sum += A(i,k)*B(k,j);
            C(i,j) = sum;
        }
    }
}

Here is the call graph for this function:

void TLinAlg::MultiplyScalar ( const double &  k,
const TFltV x,
TFltV y 
) [static]

Definition at line 414 of file linalg.cpp.

References Assert, and TVec< TVal, TSizeTy >::Len().

Referenced by Normalize(), NormalizeL1(), NormalizeLinf(), and TSparseSVD::SimpleLanczos().

                                                                      {
    Assert(x.Len() == y.Len());
    int Len = x.Len();
    for (int i = 0; i < Len; i++)
        y[i] = k * x[i];
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::MultiplyScalar ( const double &  k,
const TIntFltKdV x,
TIntFltKdV y 
) [static]

Definition at line 421 of file linalg.cpp.

References Assert, and TVec< TVal, TSizeTy >::Len().

                                                                                {
    Assert(x.Len() == y.Len());
    int Len = x.Len();
    for (int i = 0; i < Len; i++)
        y[i].Dat = k * x[i].Dat;
}

Here is the call graph for this function:

void TLinAlg::MultiplyT ( const TFltVV A,
const TFltV x,
TFltV y 
) [static]

Definition at line 468 of file linalg.cpp.

References Assert, TVVec< TVal >::GetCols(), TVVec< TVal >::GetRows(), and TVec< TVal, TSizeTy >::Len().

                                                                 {
    Assert(A.GetRows() == x.Len() && A.GetCols() == y.Len());
    int n = A.GetCols(), m = A.GetRows();
    for (int i = 0; i < n; i++) {
        y[i] = 0.0;
        for (int j = 0; j < m; j++)
            y[i] += A(j,i) * x[j];
    }
}

Here is the call graph for this function:

double TLinAlg::Norm ( const TFltV x) [static]

Definition at line 324 of file linalg.cpp.

References Norm2().

Referenced by GS(), TSparseSVD::Lanczos(), TSparseSVD::Lanczos2(), TAGMFit::MLEGradAscentGivenCAG(), TLogRegFit::MLEGradient(), Normalize(), TSparseSVD::OrtoIterSVD(), and TSparseSVD::SimpleLanczos().

                                   {
    return sqrt(Norm2(x));
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::Norm ( const TIntFltKdV x) [static]

Definition at line 341 of file linalg.cpp.

References Norm2().

                                        {
    return sqrt(Norm2(x));
}

Here is the call graph for this function:

double TLinAlg::Norm ( const TFltVV X,
int  ColId 
) [static]

Definition at line 353 of file linalg.cpp.

References Norm2().

                                               {
    return sqrt(Norm2(X, ColId));
}

Here is the call graph for this function:

double TLinAlg::Norm2 ( const TFltV x) [static]

Definition at line 320 of file linalg.cpp.

References DotProduct().

Referenced by AssertOrtogonality(), TCesna::MLEGradAscent(), TCesna::MLEGradAscentParallel(), and Norm().

                                    {
    return DotProduct(x, x);
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::Norm2 ( const TIntFltKdV x) [static]

Definition at line 333 of file linalg.cpp.

References TVec< TVal, TSizeTy >::Len(), and TMath::Sqr().

                                         {
    double Result = 0;
    for (int i = 0; i < x.Len(); i++) {
        Result += TMath::Sqr(x[i].Dat);
    }
    return Result;
}

Here is the call graph for this function:

double TLinAlg::Norm2 ( const TFltVV X,
int  ColId 
) [static]

Definition at line 349 of file linalg.cpp.

References DotProduct().

                                                {
    return DotProduct(X, ColId, X, ColId);
}

Here is the call graph for this function:

void TLinAlg::Normalize ( TFltV x) [static]

Definition at line 328 of file linalg.cpp.

References MultiplyScalar(), and Norm().

Referenced by GS().

                                {
    const double xNorm = Norm(x);
    if (xNorm > 0.0) { MultiplyScalar(1/xNorm, x, x); }
}

Here is the call graph for this function:

Here is the caller graph for this function:

void TLinAlg::Normalize ( TIntFltKdV x) [static]

Definition at line 345 of file linalg.cpp.

References MultiplyScalar(), and Norm().

                                     {
    MultiplyScalar(1/Norm(x), x, x);
}

Here is the call graph for this function:

void TLinAlg::NormalizeL1 ( TFltV x) [static]

Definition at line 380 of file linalg.cpp.

References MultiplyScalar(), and NormL1().

                                  {
    const double xNorm = NormL1(x);
    if (xNorm > 0.0) { MultiplyScalar(1/xNorm, x, x); }
}

Here is the call graph for this function:

void TLinAlg::NormalizeL1 ( TIntFltKdV x) [static]

Definition at line 385 of file linalg.cpp.

References MultiplyScalar(), and NormL1().

                                       {
    const double xNorm = NormL1(x);
    if (xNorm > 0.0) { MultiplyScalar(1/xNorm, x, x); }
}

Here is the call graph for this function:

void TLinAlg::NormalizeLinf ( TFltV x) [static]

Definition at line 404 of file linalg.cpp.

References MultiplyScalar(), and NormLinf().

                                    {
    const double xNormLinf = NormLinf(x);
    if (xNormLinf > 0.0) { MultiplyScalar(1.0/xNormLinf, x, x); }
}

Here is the call graph for this function:

void TLinAlg::NormalizeLinf ( TIntFltKdV x) [static]

Definition at line 409 of file linalg.cpp.

References MultiplyScalar(), and NormLinf().

                                         {
    const double xNormLInf = NormLinf(x);
    if (xNormLInf> 0.0) { MultiplyScalar(1.0/xNormLInf, x, x); }
}

Here is the call graph for this function:

double TLinAlg::NormL1 ( const TFltV x) [static]

Definition at line 357 of file linalg.cpp.

References TFlt::Abs(), and TVec< TVal, TSizeTy >::Len().

Referenced by NormalizeL1().

                                     {
    double norm = 0.0; const int Len = x.Len();
    for (int i = 0; i < Len; i++)
        norm += TFlt::Abs(x[i]);
    return norm;
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::NormL1 ( double  k,
const TFltV x,
const TFltV y 
) [static]

Definition at line 364 of file linalg.cpp.

References TFlt::Abs(), Assert, and TVec< TVal, TSizeTy >::Len().

                                                               {
    Assert(x.Len() == y.Len());
    double norm = 0.0; const int len = x.Len();
    for (int i = 0; i < len; i++) {
        norm += TFlt::Abs(k * x[i] + y[i]);
    }
    return norm;
}

Here is the call graph for this function:

double TLinAlg::NormL1 ( const TIntFltKdV x) [static]

Definition at line 373 of file linalg.cpp.

References TFlt::Abs(), and TVec< TVal, TSizeTy >::Len().

                                          {
    double norm = 0.0; const int Len = x.Len();
    for (int i = 0; i < Len; i++)
        norm += TFlt::Abs(x[i].Dat);
    return norm;
}

Here is the call graph for this function:

double TLinAlg::NormLinf ( const TFltV x) [static]

Definition at line 390 of file linalg.cpp.

References TFlt::Abs(), TFlt::GetMx(), and TVec< TVal, TSizeTy >::Len().

Referenced by NormalizeLinf().

                                       {
    double norm = 0.0; const int Len = x.Len();
    for (int i = 0; i < Len; i++)
        norm = TFlt::GetMx(TFlt::Abs(x[i]), norm);
    return norm;
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::NormLinf ( const TIntFltKdV x) [static]

Definition at line 397 of file linalg.cpp.

References TFlt::Abs(), TFlt::GetMx(), and TVec< TVal, TSizeTy >::Len().

                                            {
    double norm = 0.0; const int Len = x.Len();
    for (int i = 0; i < Len; i++)
        norm = TFlt::GetMx(TFlt::Abs(x[i].Dat), norm);
    return norm;
}

Here is the call graph for this function:

void TLinAlg::Rotate ( const double &  OldX,
const double &  OldY,
const double &  Angle,
double &  NewX,
double &  NewY 
) [static]

Definition at line 611 of file linalg.cpp.

                                                                                                            {
    NewX = OldX*cos(Angle) - OldY*sin(Angle);
    NewY = OldX*sin(Angle) + OldY*cos(Angle);
}
double TLinAlg::SumVec ( const TFltV x) [static]

Definition at line 279 of file linalg.cpp.

References TVec< TVal, TSizeTy >::Len().

Referenced by TAGMFit::Likelihood().

                                     {
    const int len = x.Len();
    double Res = 0.0;
    for (int i = 0; i < len; i++) {
        Res += x[i];
    }
    return Res;
}

Here is the call graph for this function:

Here is the caller graph for this function:

double TLinAlg::SumVec ( double  k,
const TFltV x,
const TFltV y 
) [static]

Definition at line 288 of file linalg.cpp.

References Assert, and TVec< TVal, TSizeTy >::Len().

                                                               {
    Assert(x.Len() == y.Len());
    const int len = x.Len();
    double Res = 0.0;
    for (int i = 0; i < len; i++) {
        Res += k * x[i] + y[i];
    }
    return Res;
}

Here is the call graph for this function:

void TLinAlg::Transpose ( const TFltVV A,
TFltVV B 
) [static]

Definition at line 539 of file linalg.cpp.

References Assert, TVVec< TVal >::At(), TVVec< TVal >::GetCols(), and TVVec< TVal >::GetRows().

                                                  {
        Assert(B.GetRows() == A.GetCols() && B.GetCols() == A.GetRows());
        for (int i = 0; i < A.GetCols(); i++) {
                for (int j = 0; j < A.GetRows(); j++) {
                        B.At(i, j) = A.At(j, i);
                }
        }
}

Here is the call graph for this function:


The documentation for this class was generated from the following files: