SNAP Library 6.0, User Reference  2020-12-09 16:24:20
SNAP, a general purpose, high performance system for analysis and manipulation of large networks
TCorr Class Reference

#include <xmath.h>

Public Member Functions

 TCorr ()
 
 TCorr (const TFltV &ValV1, const TFltV &ValV2)
 
 TCorr (TSIn &)
 
void Save (TSOut &)
 
TCorroperator= (const TCorr &)
 
double GetCorrCf () const
 
double GetCorrCfPrb () const
 
TStr GetStr () const
 

Static Public Member Functions

static PCorr New (const TFltV &ValV1, const TFltV &ValV2)
 
static PCorr Load (TSIn &SIn)
 

Private Attributes

TCRef CRef
 
int ValVLen
 
double CorrCf
 
double CorrCfPrb
 
double FisherZ
 

Friends

class TPt< TCorr >
 

Detailed Description

Definition at line 269 of file xmath.h.

Constructor & Destructor Documentation

TCorr::TCorr ( )
inline

Definition at line 276 of file xmath.h.

276 {}
TCorr::TCorr ( const TFltV ValV1,
const TFltV ValV2 
)

Definition at line 554 of file xmath.cpp.

554  :
555  ValVLen(ValV1.Len()), CorrCf(), CorrCfPrb(), FisherZ(){
556  static const double TINY=1.0e-20;
557  IAssert(ValV1.Len()==ValV2.Len());
558 
559  // calculate the means
560  double MeanVal1=0; double MeanVal2=0;
561  {for (int ValN=0; ValN<ValVLen; ValN++){
562  MeanVal1+=ValV1[ValN]; MeanVal2+=ValV2[ValN];}}
563  MeanVal1/=ValVLen; MeanVal2/=ValVLen;
564 
565  // calculate correlation coefficient
566  double yt, xt;
567  double syy=0.0; double sxy=0.0; double sxx=0.0;
568  {for (int ValN=0; ValN<ValVLen; ValN++){
569  xt=ValV1[ValN]-MeanVal1;
570  yt=ValV2[ValN]-MeanVal2;
571  sxx+=xt*xt;
572  syy+=yt*yt;
573  sxy+=xt*yt;
574  }}
575  if (sxx*syy==0){
576  CorrCf=0; //** not in numerical recipes - check why (pojavi se, ko so same nicle)
577  } else {
578  CorrCf=sxy/sqrt(sxx*syy);
579  }
580  // calculate correlation coefficient significance level
581  double df=ValVLen-2;
582  double t=CorrCf*sqrt(df/((1.0-CorrCf+TINY)*(1.0+CorrCf+TINY)));
583  CorrCfPrb=TSpecFunc::BetaI(0.5*df,0.5,df/(df+t*t));
584  // calculate Fisher's Z transformation
585  FisherZ=0.5*log((1.0+(CorrCf)+TINY)/(1.0-(CorrCf)+TINY));
586 }
#define IAssert(Cond)
Definition: bd.h:262
double FisherZ
Definition: xmath.h:274
double CorrCf
Definition: xmath.h:272
TSizeTy Len() const
Returns the number of elements in the vector.
Definition: ds.h:575
int ValVLen
Definition: xmath.h:271
static double BetaI(const double &a, const double &b, const double &x)
Definition: xmath.cpp:137
double CorrCfPrb
Definition: xmath.h:273
TCorr::TCorr ( TSIn )
inline

Definition at line 280 of file xmath.h.

280 {Fail;}
#define Fail
Definition: bd.h:238

Member Function Documentation

double TCorr::GetCorrCf ( ) const
inline

Definition at line 286 of file xmath.h.

286 {return CorrCf;}
double CorrCf
Definition: xmath.h:272
double TCorr::GetCorrCfPrb ( ) const
inline

Definition at line 287 of file xmath.h.

287 {return CorrCfPrb;}
double CorrCfPrb
Definition: xmath.h:273
TStr TCorr::GetStr ( ) const
static PCorr TCorr::Load ( TSIn SIn)
inlinestatic

Definition at line 281 of file xmath.h.

281 {return new TCorr(SIn);}
TCorr()
Definition: xmath.h:276
static PCorr TCorr::New ( const TFltV ValV1,
const TFltV ValV2 
)
inlinestatic

Definition at line 278 of file xmath.h.

278  {
279  return PCorr(new TCorr(ValV1, ValV2));}
TCorr()
Definition: xmath.h:276
TPt< TCorr > PCorr
Definition: xmath.h:269
TCorr& TCorr::operator= ( const TCorr )
inline

Definition at line 284 of file xmath.h.

284 {Fail; return *this;}
#define Fail
Definition: bd.h:238
void TCorr::Save ( TSOut )
inline

Definition at line 282 of file xmath.h.

282 {Fail;}
#define Fail
Definition: bd.h:238

Friends And Related Function Documentation

friend class TPt< TCorr >
friend

Definition at line 269 of file xmath.h.

Member Data Documentation

double TCorr::CorrCf
private

Definition at line 272 of file xmath.h.

double TCorr::CorrCfPrb
private

Definition at line 273 of file xmath.h.

TCRef TCorr::CRef
private

Definition at line 269 of file xmath.h.

double TCorr::FisherZ
private

Definition at line 274 of file xmath.h.

int TCorr::ValVLen
private

Definition at line 271 of file xmath.h.


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