| 
    SNAP Library 2.0, Developer Reference
    2013-05-13 16:33:57
    
   SNAP, a general purpose, high performance system for analysis and manipulation of large networks 
   | 
  
  
  
 
#include <dt.h>
Public Member Functions | |
| TRnd (const int &_Seed=1, const int &Steps=0) | |
| TRnd (TSIn &SIn) | |
| void | Save (TSOut &SOut) const | 
| void | LoadXml (const PXmlTok &XmlTok, const TStr &Nm) | 
| void | SaveXml (TSOut &SOut, const TStr &Nm) const | 
| TRnd & | operator= (const TRnd &Rnd) | 
| bool | operator== (const TRnd &) const | 
| double | GetUniDev () | 
| int | GetUniDevInt (const int &Range=0) | 
| int | GetUniDevInt (const int &MnVal, const int &MxVal) | 
| uint | GetUniDevUInt (const uint &Range=0) | 
| int64 | GetUniDevInt64 (const int64 &Range=0) | 
| uint64 | GetUniDevUInt64 (const uint64 &Range=0) | 
| double | GetNrmDev () | 
| double | GetNrmDev (const double &Mean, const double &SDev, const double &Mn, const double &Mx) | 
| double | GetExpDev () | 
| double | GetExpDev (const double &Lambda) | 
| double | GetGammaDev (const int &Order) | 
| double | GetPoissonDev (const double &Mean) | 
| double | GetBinomialDev (const double &Prb, const int &Trials) | 
| int | GetGeoDev (const double &Prb) | 
| double | GetPowerDev (const double &AlphaSlope) | 
| double | GetRayleigh (const double &Sigma) | 
| double | GetWeibull (const double &K, const double &Lambda) | 
| void | PutSeed (const int &_Seed) | 
| int | GetSeed () const | 
| void | Randomize () | 
| void | Move (const int &Steps) | 
| bool | Check () | 
| void | SaveTxt (TOLx &Lx) const | 
Static Public Member Functions | |
| static double | GetUniDevStep (const int &Seed, const int &Steps) | 
| static double | GetNrmDevStep (const int &Seed, const int &Steps) | 
| static double | GetExpDevStep (const int &Seed, const int &Steps) | 
| static TRnd | LoadTxt (TILx &Lx) | 
Static Public Attributes | |
| static const int | RndSeed = 0 | 
Private Member Functions | |
| int | GetNextSeed () | 
Private Attributes | |
| int | Seed | 
Static Private Attributes | |
| static const int | a = 16807 | 
| static const int | m = 2147483647 | 
| static const int | q = 127773 | 
| static const int | r = 2836 | 
| TRnd::TRnd | ( | const int & | _Seed = 1,  | 
        
| const int & | Steps = 0  | 
        ||
| ) |  [inline] | 
        
| TRnd::TRnd | ( | TSIn & | SIn | ) |  [inline, explicit] | 
        
| bool TRnd::Check | ( | ) | 
Definition at line 33 of file dt.cpp.
References GetNextSeed(), and Seed.
                {
  int PSeed=Seed; Seed=1;
  for (int SeedN=0; SeedN<10000; SeedN++){GetNextSeed();}
  bool Ok=Seed==1043618065; Seed=PSeed; return Ok;
}

| double TRnd::GetBinomialDev | ( | const double & | Prb, | 
| const int & | Trials | ||
| ) | 
Definition at line 154 of file dt.cpp.
References GetUniDev(), TSpecFunc::LnGamma(), and TMath::Pi.
Referenced by TAGM::RndConnectInsideCommunity().
                                                               {
  int j;
  static int nold=(-1);
  double am,em,g,angle,p,bnl,sq,t,y;
  static double pold=(-1.0),pc,plog,pclog,en,oldg;
  p=(Prb <= 0.5 ? Prb : 1.0-Prb);
  am=Trials*p;
  if (Trials < 25) {
    bnl=0.0;
    for (j=1;j<=Trials;j++)
      if (GetUniDev() < p) ++bnl;
  } else if (am < 1.0) {
    g=exp(-am);
    t=1.0;
    for (j=0;j<=Trials;j++) {
      t *= GetUniDev();
      if (t < g) break;
    }
    bnl=(j <= Trials ? j : Trials);
  } else {
    if (Trials != nold) {
      en=Trials;
      oldg=TSpecFunc::LnGamma(en+1.0);
      nold=Trials;
    } if (p != pold) {
      pc=1.0-p;
      plog=log(p);
      pclog=log(pc);
      pold=p;
    }
    sq=sqrt(2.0*am*pc);
    do {
      do {
        angle=TMath::Pi*GetUniDev();
        y=tan(angle);
        em=sq*y+am;
      } while (em < 0.0 || em >= (en+1.0));
      em=floor(em);
      t=1.2*sq*(1.0+y*y)*exp(oldg-(em+1.0)
        -TSpecFunc::LnGamma(en-em+1.0)+em*plog+(en-em)*pclog);
    } while (GetUniDev() > t);
    bnl=em;
  }
  if (p != Prb) bnl=Trials-bnl;
  return bnl;
}


| double TRnd::GetExpDev | ( | ) | 
Definition at line 83 of file dt.cpp.
References GetUniDev().
Referenced by TNetInfBs::GenCascade(), TNIBs::GenCascade(), GetExpDev(), and GetExpDevStep().
                      {
  double UniDev;
  do {
    UniDev=GetUniDev();
  } while (UniDev==0.0);
  return -log(UniDev);
}


| double TRnd::GetExpDev | ( | const double & | Lambda | ) | 
Definition at line 91 of file dt.cpp.
References GetExpDev().
                                           {
  return GetExpDev()/Lambda;
}

| static double TRnd::GetExpDevStep | ( | const int & | Seed, | 
| const int & | Steps | ||
| ) |  [inline, static] | 
        
Definition at line 68 of file dt.h.
References GetExpDev(), and Move().
Referenced by TExpBi::GetBiFuncVal().


| double TRnd::GetGammaDev | ( | const int & | Order | ) | 
Definition at line 95 of file dt.cpp.
References Fail, and GetUniDev().
Referenced by TMAGNodeBeta::AttrGen().
                                        {
  int j;
  double am,e,s,v1,v2,x,y;
  if (Order<1){Fail;}
  if (Order<6) {
    x=1.0;
    for (j=1;j<=Order;j++) x *=GetUniDev();
    x = -log(x);
  } else {
    do {
      do {
        do {
          v1=2.0*GetUniDev()-1.0;
          v2=2.0*GetUniDev()-1.0;
        } while (v1*v1+v2*v2 > 1.0);
        y=v2/v1;
        am=Order-1;
        s=sqrt(2.0*am+1.0);
        x=s*y+am;
      } while (x <= 0.0);
      e=(1.0+y*y)*exp(am*log(x/am)-s*y);
    } while (GetUniDev()>e);
  }
  return x;
}


| int TRnd::GetGeoDev | ( | const double & | Prb | ) |  [inline] | 
        
Definition at line 45 of file dt.h.
References GetUniDev().
Referenced by TForestFire::BurnGeoFire(), and TUndirFFire::BurnGeoFire().
                                  {
    return 1+(int)floor(log(1.0-GetUniDev())/log(1.0-Prb));}


| int TRnd::GetNextSeed | ( | ) |  [inline, private] | 
        
| double TRnd::GetNrmDev | ( | ) | 
Definition at line 63 of file dt.cpp.
References GetUniDev().
Referenced by TMAGAffMtx::AddRndNoise(), TKronMtx::AddRndNoise(), TLAMisc::FillRnd(), GetNrmDev(), GetNrmDevStep(), and TSnap::TSnapDetail::GetSphereDev().
                      {
  double v1, v2, rsq;
  do {
    v1=2.0*GetUniDev()-1.0; // pick two uniform numbers in the square
    v2=2.0*GetUniDev()-1.0; // extending from -1 to +1 in each direction
    rsq=v1*v1+v2*v2; // see if they are in the unit cicrcle
  } while ((rsq>=1.0)||(rsq==0.0)); // and if they are not, try again
  double fac=sqrt(-2.0*log(rsq)/rsq); // Box-Muller transformation
  return v1*fac;
//  return v2*fac; // second deviate
}


| double TRnd::GetNrmDev | ( | const double & | Mean, | 
| const double & | SDev, | ||
| const double & | Mn, | ||
| const double & | Mx | ||
| ) | 
Definition at line 75 of file dt.cpp.
References GetNrmDev().
                                                                            {
  double Val=Mean+GetNrmDev()*SDev;
  if (Val<Mn){Val=Mn;}
  if (Val>Mx){Val=Mx;}
  return Val;
}

| static double TRnd::GetNrmDevStep | ( | const int & | Seed, | 
| const int & | Steps | ||
| ) |  [inline, static] | 
        
Definition at line 66 of file dt.h.
References GetNrmDev(), and Move().
Referenced by TExpBi::GetBiFuncVal().


| double TRnd::GetPoissonDev | ( | const double & | Mean | ) | 
Definition at line 121 of file dt.cpp.
References GetUniDev(), TSpecFunc::LnGamma(), and TMath::Pi.
                                            {
  static double sq,alxm,g,oldm=(-1.0);
  double em,t,y;
  if (Mean < 12.0) {
    if (Mean != oldm) {
      oldm=Mean;
      g=exp(-Mean);
    }
    em = -1;
    t=1.0;
    do {
      ++em;
      t *= GetUniDev();
    } while (t>g);
  } else {
    if (Mean != oldm) {
      oldm=Mean;
      sq=sqrt(2.0*Mean);
      alxm=log(Mean);
      g=Mean*alxm-TSpecFunc::LnGamma(Mean+1.0);
    }
    do {
      do {
        y=tan(TMath::Pi*GetUniDev());
        em=sq*y+Mean;
      } while (em < 0.0);
      em=floor(em);
      t=0.9*(1.0+y*y)*exp(em*alxm-TSpecFunc::LnGamma(em+1.0)-g);
    } while (GetUniDev()>t);
  }
  return em;
}

| double TRnd::GetPowerDev | ( | const double & | AlphaSlope | ) |  [inline] | 
        
Definition at line 47 of file dt.h.
References GetUniDev(), and IAssert.
Referenced by TNetInfBs::GenCascade(), TNIBs::GenCascade(), TAGMUtil::GenPLSeq(), and TSnap::GenRndPowerLaw().
                                              { // power-law degree distribution (AlphaSlope>0)
    IAssert(AlphaSlope>1.0);
    return pow(1.0-GetUniDev(), -1.0/(AlphaSlope-1.0));}


| double TRnd::GetRayleigh | ( | const double & | Sigma | ) |  [inline] | 
        
Definition at line 50 of file dt.h.
References GetUniDev(), and IAssert.
Referenced by TNetInfBs::GenCascade(), and TNIBs::GenCascade().


| int TRnd::GetSeed | ( | ) |  const [inline] | 
        
| double TRnd::GetUniDev | ( | ) |  [inline] | 
        
Definition at line 30 of file dt.h.
References GetNextSeed(), and m.
Referenced by TFfGGen::AddNodes(), TForestFire::BurnExpFire(), TMAGFitBern::ComputeJointLL(), TMAGFitBern::DoEMAlg(), TNIBs::GenCascade(), TSnap::GenCopyModel(), TKronMtx::GenDetKronecker(), TKronMtx::GenFastKronecker(), TKronMtx::GenKronecker(), TSnap::GenRMat(), TKronMtx::GenRndGraph(), TSnap::GenSmallWorld(), GetBinomialDev(), GetExpDev(), GetGammaDev(), GetGeoDev(), TKronMtx::GetInitMtx(), GetNrmDev(), GetPoissonDev(), GetPowerDev(), GetRayleigh(), TFlt::GetRnd(), GetUniDevStep(), GetWeibull(), TCluster::MCMC(), TAGMFast::NeighborComInit(), TSparseSVD::OrtoIterSVD(), TAGMFast::RandomInit(), TMAGFitBern::RandomInit(), TAGMFit::RunMCMC(), TKroneckerLL::SampleNextPerm(), TKroneckerLL::SetRandomEdges(), TMAGAffMtx::SetRndMtx(), and TKronMtx::SetRndMtx().
{return GetNextSeed()/double(m);}


| int TRnd::GetUniDevInt | ( | const int & | Range = 0 | ) | 
Definition at line 39 of file dt.cpp.
References GetNextSeed(), and Seed.
Referenced by TFfGGen::AddNodes(), TUndirFFire::AddNodes(), TNIBs::BSG(), TAGMUtil::ConnectCmtyVV(), TMAGFitBern::DoEStepApxOneIter(), TMAGFitBern::DoEStepOneIter(), TAGMFast::FindComsByCV(), TKronNoise::FlipEdgeNoise(), TSnap::GenConfModel(), TSnap::GenCopyModel(), TSnap::GenGeoPrefAttach(), TTree< TVal >::GenRandomTree(), TSnap::GenRndBipart(), TSnap::GenRndGnm(), TSnap::GenSmallWorld(), TTimeNENet::GetGnmRndNet(), TSnap::GetMxDegNId(), TSnap::GetMxInDegNId(), TSnap::GetMxOutDegNId(), TTimeNENet::GetPrefAttach(), TBool::GetRnd(), TInt::GetRnd(), TSnap::TSnapDetail::GetRndEdgeNonAdjNode(), THash< TKey, TDat, THashFunc >::GetRndKeyId(), THashSet< TInt >::GetRndKeyId(), TCnCom::GetRndNId(), TBPGraph::GetRndNId(), TBreathFS< PGraph >::GetRndPath(), GetUniDevInt(), GetUniDevInt64(), GetUniDevUInt64(), TKroneckerLL::MetroGibbsSampleNext(), TAGMFast::NeighborComInit(), TAGMFast::RandomInit(), TAGMFit::RandomInit(), TAGMUtil::RewireCmtyNID(), TAGM::RndConnectInsideCommunity(), TKroneckerLL::SampleNextPerm(), TAGMFit::SampleTransition(), TNIBs::SG(), TVec< TVal, TSizeTy >::Shuffle(), TGLib_OLD::TVec< TVal >::Shuffle(), TVecPool< TVal, TSizeTy >::ShuffleAll(), TGLib_OLD::TVecPool< TVal >::ShuffleAll(), TVVec< TVal >::ShuffleX(), TVVec< TVal >::ShuffleY(), TUniChDb::TestFindNextWordOrSentenceBoundary(), and TCluster::Train().
                                      {
  int Seed=GetNextSeed();
  if (Range==0){return Seed;}
  else {return Seed%Range;}
}


| int TRnd::GetUniDevInt | ( | const int & | MnVal, | 
| const int & | MxVal | ||
| ) |  [inline] | 
        
Definition at line 32 of file dt.h.
References GetUniDevInt(), and IAssert.
                                                      {
    IAssert(MnVal<=MxVal); return MnVal+GetUniDevInt(MxVal-MnVal+1);}

| int64 TRnd::GetUniDevInt64 | ( | const int64 & | Range = 0 | ) | 
Definition at line 51 of file dt.cpp.
References GetUniDevInt().
Referenced by TVec< TVal, TSizeTy >::Shuffle().
                                            {
  const int64 RndVal = int64((uint64(GetUniDevInt())<<32) | uint64(GetUniDevInt()));
  if (Range==0){return RndVal;}
  else {return RndVal%Range;}
}


| static double TRnd::GetUniDevStep | ( | const int & | Seed, | 
| const int & | Steps | ||
| ) |  [inline, static] | 
        
Definition at line 64 of file dt.h.
References GetUniDev(), and Move().
Referenced by TExpBi::GetBiFuncVal().


| uint TRnd::GetUniDevUInt | ( | const uint & | Range = 0 | ) | 
Definition at line 45 of file dt.cpp.
References GetNextSeed(), and Seed.
Referenced by TUInt::GetRnd(), and TUniCodec::GetRndUint().
                                         {
  uint Seed=uint(GetNextSeed()%0x10000)*0x10000+uint(GetNextSeed()%0x10000);
  if (Range==0){return Seed;}
  else {return Seed%Range;}
}


| uint64 TRnd::GetUniDevUInt64 | ( | const uint64 & | Range = 0 | ) | 
Definition at line 57 of file dt.cpp.
References GetUniDevInt().
                                               {
 const uint64 RndVal = uint64((uint64(GetUniDevInt())<<32) | uint64(GetUniDevInt()));
 if (Range==0){return RndVal;}
 else {return RndVal%Range;}
}

| double TRnd::GetWeibull | ( | const double & | K, | 
| const double & | Lambda | ||
| ) |  [inline] | 
        
| TRnd TRnd::LoadTxt | ( | TILx & | Lx | ) |  [static] | 
        
| void TRnd::LoadXml | ( | const PXmlTok & | XmlTok, | 
| const TStr & | Nm | ||
| ) | 
Definition at line 9 of file dt.cpp.
References TXmlObjSer::GetIntArg(), Seed, and XLoadHd.
                                                       {
  XLoadHd(Nm);
  Seed=TXmlObjSer::GetIntArg(XmlTok, "Seed");
}

| void TRnd::Move | ( | const int & | Steps | ) | 
Definition at line 29 of file dt.cpp.
References GetNextSeed().
Referenced by GetExpDevStep(), GetNrmDevStep(), GetUniDevStep(), and TRnd().
                               {
  for (int StepN=0; StepN<Steps; StepN++){GetNextSeed();}
}


| bool TRnd::operator== | ( | const TRnd & | ) |  const [inline] | 
        
| void TRnd::PutSeed | ( | const int & | _Seed | ) | 
Definition at line 18 of file dt.cpp.
Referenced by TAGMFit::Load(), TAGMFast::Load(), TKronMtx::PutRndSeed(), TMAGFitBern::RandomInit(), Randomize(), and TRnd().
                                  {
  Assert(_Seed>=0);
  if (_Seed==0){
    //Seed=int(time(NULL));
    Seed=abs(int(TSysTm::GetPerfTimerTicks()));
  } else {
    Seed=_Seed;
    //Seed=abs(_Seed*100000)+1;
  }
}

| void TRnd::Randomize | ( | ) |  [inline] | 
        
| void TRnd::Save | ( | TSOut & | SOut | ) |  const [inline] | 
        
| void TRnd::SaveTxt | ( | TOLx & | Lx | ) | const | 
| void TRnd::SaveXml | ( | TSOut & | SOut, | 
| const TStr & | Nm | ||
| ) | const | 
Definition at line 14 of file dt.cpp.
References TInt::GetStr(), Seed, and XSaveBETagArg.
                                                    {
  XSaveBETagArg(Nm, "Seed", TInt::GetStr(Seed));
}

const int TRnd::a = 16807 [static, private] | 
        
Definition at line 15 of file dt.h.
Referenced by GetNextSeed().
const int TRnd::m = 2147483647 [static, private] | 
        
Definition at line 15 of file dt.h.
Referenced by GetNextSeed(), and GetUniDev().
const int TRnd::q = 127773 [static, private] | 
        
Definition at line 15 of file dt.h.
Referenced by GetNextSeed().
const int TRnd::r = 2836 [static, private] | 
        
Definition at line 15 of file dt.h.
Referenced by GetNextSeed().
const int TRnd::RndSeed = 0 [static] | 
        
Definition at line 13 of file dt.h.
Referenced by Randomize().
int TRnd::Seed [private] | 
        
Definition at line 16 of file dt.h.
Referenced by Check(), GetNextSeed(), GetSeed(), GetUniDevInt(), GetUniDevUInt(), LoadXml(), operator=(), PutSeed(), Save(), SaveTxt(), SaveXml(), and TRnd().