SNAP Library, Developer Reference  2012-10-02 12:56:23
SNAP, a general purpose network analysis and graph mining library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
xmath.cpp
Go to the documentation of this file.
00001 
00002 // Mathematical-Utilities
00003 double TMath::E=2.71828182845904523536;
00004 double TMath::Pi=3.14159265358979323846;
00005 double TMath::LogOf2=log(double(2));
00006 
00008 // Special-Functions
00009 void TSpecFunc::GammaPSeries/*gser*/(
00010  double& gamser, const double& a, const double& x, double& gln){
00011   static const int ITMAX=100;
00012   static const double EPS=3.0e-7;
00013   int n;
00014   double sum, del, ap;
00015 
00016   gln=LnGamma(a);
00017   if (x <= 0.0){
00018     IAssert(x>=0); /*if (x < 0.0) nrerror("x less than 0 in routine gser");*/
00019     gamser=0.0;
00020     return;
00021   } else {
00022     ap=a;
00023     del=sum=1.0/a;
00024     for (n=1; n<=ITMAX; n++){
00025       ++ap;
00026       del *= x/ap;
00027       sum += del;
00028       if (fabs(del) < fabs(sum)*EPS){
00029         gamser=sum*exp(-x+a*log(x)-(gln));
00030         return;
00031       }
00032     }
00033     Fail; /*nrerror("a too large, ITMAX too small in routine gser");*/
00034     return;
00035   }
00036 }
00037 
00038 void TSpecFunc::GammaQContFrac/*gcf*/(
00039  double& gammcf, const double& a, const double& x, double& gln){
00040   static const int ITMAX=100;
00041   static const double EPS=3.0e-7;
00042   static const double  FPMIN=1.0e-30;
00043   int i;
00044   double an, b, c, d, del, h;
00045 
00046   gln=LnGamma(a);
00047   b=x+1.0-a;
00048   c=1.0/FPMIN;
00049   d=1.0/b;
00050   h=d;
00051   for (i=1;i<=ITMAX;i++){
00052     an = -i*(i-a);
00053     b += 2.0;
00054     d=an*d+b;
00055     if (fabs(d) < FPMIN) d=FPMIN;
00056     c=b+an/c;
00057     if (fabs(c) < FPMIN) c=FPMIN;
00058     d=1.0/d;
00059     del=d*c;
00060     h *= del;
00061     if (fabs(del-1.0) < EPS) break;
00062   }
00063   IAssert(i<=ITMAX);
00064   /*if (i > ITMAX) nrerror("a too large, ITMAX too small in gcf");*/
00065   gammcf=exp(-x+a*log(x)-(gln))*h;
00066 }
00067 
00068 double TSpecFunc::GammaQ/*gammq*/(const double& a, const double& x){
00069   IAssert((x>=0)&&(a>0));
00070   double gamser, gammcf, gln;
00071   if (x<(a+1.0)){
00072     GammaPSeries(gamser,a,x,gln);
00073     return 1.0-gamser;
00074   } else {
00075     GammaQContFrac(gammcf,a,x,gln);
00076     return gammcf;
00077   }
00078 }
00079 
00080 double TSpecFunc::LnGamma/*gammln*/(const double& xx){
00081   double x, y, tmp, ser;
00082   static double cof[6]={76.18009172947146,-86.50532032941677,
00083           24.01409824083091,-1.231739572450155,
00084           0.1208650973866179e-2,-0.5395239384953e-5};
00085   int j;
00086 
00087   y=x=xx;
00088   tmp=x+5.5;
00089   tmp -= (x+0.5)*log(tmp);
00090   ser=1.000000000190015;
00091   for (j=0;j<=5;j++) ser += cof[j]/++y;
00092   return -tmp+log(2.5066282746310005*ser/x);
00093 }
00094 
00095 double TSpecFunc::LnComb(const int& n, const int& k){
00096   return LnGamma(n+1)-LnGamma(k+1)-LnGamma(n-k+1);
00097 }
00098 
00099 double TSpecFunc::BetaCf(const double& a, const double& b, const double& x){
00100   static const double MAXIT=100;
00101   static const double EPS=3.0e-7;
00102   static const double FPMIN=1.0e-30;
00103   int m,m2;
00104   double aa,c,d,del,h,qab,qam,qap;
00105 
00106   qab=a+b;
00107   qap=a+1.0;
00108   qam=a-1.0;
00109   c=1.0;
00110   d=1.0-qab*x/qap;
00111   if (fabs(d) < FPMIN) d=FPMIN;
00112   d=1.0/d;
00113   h=d;
00114   for (m=1;m<=MAXIT;m++) {
00115     m2=2*m;
00116     aa=m*(b-m)*x/((qam+m2)*(a+m2));
00117     d=1.0+aa*d;
00118     if (fabs(d) < FPMIN) d=FPMIN;
00119     c=1.0+aa/c;
00120     if (fabs(c) < FPMIN) c=FPMIN;
00121     d=1.0/d;
00122     h *= d*c;
00123     aa = -(a+m)*(qab+m)*x/((a+m2)*(qap+m2));
00124     d=1.0+aa*d;
00125     if (fabs(d) < FPMIN) d=FPMIN;
00126     c=1.0+aa/c;
00127     if (fabs(c) < FPMIN) c=FPMIN;
00128     d=1.0/d;
00129     del=d*c;
00130     h *= del;
00131     if (fabs(del-1.0) < EPS) break;
00132   }
00133   if (m > MAXIT){Fail;}// a or b too big, or MAXIT too small in betacf
00134   return h;
00135 }
00136 
00137 double TSpecFunc::BetaI(const double& a, const double& b, const double& x){
00138   double bt;
00139 
00140   if (x < 0.0 || x > 1.0){Fail;} // Bad x in routine betai
00141   if (x == 0.0 || x == 1.0) bt=0.0;
00142   else
00143     bt=exp(LnGamma(a+b)-LnGamma(a)-LnGamma(b)+a*log(x)+b*log(1.0-x));
00144   if (x < (a+1.0)/(a+b+2.0))
00145     return bt*BetaCf(a,b,x)/a;
00146   else
00147     return 1.0-bt*BetaCf(b,a,1.0-x)/b;
00148 }
00149 
00150 void TSpecFunc::LinearFit(
00151  const TVec<TFltPr>& XY, double& A, double& B,
00152  double& SigA, double& SigB, double& Chi2, double& R2) {
00153   // y = a + bx :: SigA (SigB): A's uncertainty; Chi2: std dev on all points
00154   int i;
00155   double t, sxoss, sx = 0.0, sy = 0.0, st2 = 0.0, ss, sigdat;
00156 
00157   A = B = SigA = SigB = Chi2 = 0.0;
00158   for (i = 0; i < XY.Len(); i++) {
00159     sx += XY[i].Val1;
00160     sy += XY[i].Val2;
00161   }
00162   ss = XY.Len();
00163   sxoss = sx / ss;
00164   for (i = 0; i <XY.Len(); i++) {
00165     t = XY[i].Val1 - sxoss;
00166     st2 += t*t;
00167     B += t * XY[i].Val2;
00168   }
00169   B /= st2;
00170   A = (sy - sx * B) / ss;
00171   SigA = sqrt((1.0 + sx * sx / (ss * st2)) / ss);
00172   SigB = sqrt(1.0 / st2);
00173   for (i = 0; i < XY.Len(); i++)
00174     Chi2 += TMath::Sqr(XY[i].Val2 - A - B * XY[i].Val1);
00175   sigdat = sqrt(Chi2 / (XY.Len() - 2));
00176   SigA *= sigdat;
00177   SigB *= sigdat;
00178 
00179   // calculate R squared
00180   { double N = XY.Len(), sXY=0.0, sX=0.0, sY=0.0, sSqX=0.0, sSqY=0.0;
00181   for (int s =0; s < XY.Len(); s++) {
00182     sX += XY[s].Val1;  sY += XY[s].Val2;
00183     sXY += XY[s].Val1 * XY[s].Val2;
00184     sSqX += TMath::Sqr(XY[s].Val1);
00185     sSqY += TMath::Sqr(XY[s].Val2);
00186   }
00187   R2 = TMath::Sqr(N*sXY - sX*sY) / ((N*sSqX - sX*sX) * (N*sSqY - sY*sY)); }
00188   if (1.1 < R2 || -1.1 > R2) R2 = 0.0;
00189   if (_isnan(A) || ! _finite(A)) A = 0.0;
00190   if (_isnan(B) || ! _finite(B)) B = 0.0;
00191 }
00192 
00193 void TSpecFunc::PowerFit(const TVec<TFltPr>& XY, double& A, double& B,
00194  double& SigA, double& SigB, double& Chi2, double& R2) {
00195   // y = a x^b :: SigA (SigB): A's uncertainty; Chi2: std dev on all points
00196   // log fit
00197   double AA, BB;
00198   TFltPrV LogXY(XY.Len(), 0);
00199   for (int s = 0; s < XY.Len(); s++) {
00200     LogXY.Add(TFltPr(log((double)XY[s].Val1), log((double)XY[s].Val2)));
00201   }
00202   TSpecFunc::LinearFit(LogXY, AA, BB, SigA, SigB, Chi2, R2);
00203   A = exp(AA);  B = BB;
00204   if (_isnan(AA) || ! _finite(AA)) A = 0.0;
00205   if (_isnan(BB) || ! _finite(BB)) B = 0.0;
00206 }
00207 
00208 void TSpecFunc::LogFit(const TVec<TFltPr>& XY, double& A, double& B,
00209  double& SigA, double& SigB, double& Chi2, double& R2) {
00210   // Y = A + B*log(X)
00211   TFltPrV LogXY(XY.Len(), 0);
00212   for (int s = 0; s < XY.Len(); s++) {
00213     LogXY.Add(TFltPr(log((double)XY[s].Val1), XY[s].Val2));
00214   }
00215   TSpecFunc::LinearFit(LogXY, A, B, SigA, SigB, Chi2, R2);
00216 }
00217 
00218 void TSpecFunc::ExpFit(const TVec<TFltPr>& XY, double& A, double& B,
00219  double& SigA, double& SigB, double& Chi2, double& R2) {
00220   // Y = A * exp(B*X)
00221   TFltPrV XLogY(XY.Len(), 0);
00222   double AA, BB;
00223   for (int s = 0; s < XY.Len(); s++) {
00224     XLogY.Add(TFltPr(XY[s].Val1, log((double)XY[s].Val2)));
00225   }
00226   TSpecFunc::LinearFit(XLogY, AA, BB, SigA, SigB, Chi2, R2);
00227   A = exp(AA);
00228   B = BB;
00229 }
00230 
00231 double TSpecFunc::Entropy(const TIntV& ValV) {
00232   TFltV NewValV(ValV.Len());
00233   for (int i = 0; i < ValV.Len(); i++) { NewValV[i] = ValV[i]; }
00234   return Entropy(NewValV);
00235 }
00236 
00237 // Entropy of a distribution: ValV[i] contains the number of events i
00238 double TSpecFunc::Entropy(const TFltV& ValV) {
00239   double Sum=0, Ent=0;
00240   for (int i = 0; i < ValV.Len(); i++) {
00241     const double& Val = ValV[i];
00242     if (Val > 0.0) { Ent -= Val * log(Val);  Sum += Val; }
00243   }
00244   if (Sum > 0.0) {
00245     Ent /= Sum;
00246     Ent += log(Sum);
00247     Ent /= TMath::LogOf2;
00248   } else return 1.0;
00249   return Ent;
00250 }
00251 
00252 void TSpecFunc::EntropyFracDim(const TIntV& ValV, TFltV& EntropyV) {
00253   TFltV NewValV(ValV.Len());
00254   for (int i = 0; i < ValV.Len(); i++) { 
00255     IAssert(ValV[i]==1 || ValV[i] == 0);
00256     NewValV[i] = ValV[i]; 
00257   }
00258   EntropyFracDim(NewValV, EntropyV);
00259 }
00260 
00261 // Entropy fractal dimension. Input is a vector {0,1}^n. 
00262 // Where 0 means the event did not occur, and 1 means it occured.
00263 // Works exactly as Mengzi Wang's code.
00264 void TSpecFunc::EntropyFracDim(const TFltV& ValV, TFltV& EntropyV) {
00265   TFltV ValV1, ValV2;
00266   int Pow2 = 1;
00267   while (2*Pow2 <= ValV.Len()) { Pow2 *= 2; }
00268   ValV1.Gen(Pow2);
00269   for (int i = 0; i < Pow2; i++) { ValV1[i] = ValV[i]; 
00270     IAssert(ValV[i]==1.0 || ValV[i] == 0.0); }
00271   EntropyV.Clr();
00272   EntropyV.Add(Entropy(ValV1)); // 2^Pow2 windows
00273   while (ValV1.Len() > 2) {
00274     ValV2.Gen(ValV1.Len() / 2);
00275     for (int i = 0; i < ValV1.Len(); i++) {
00276       ValV2[i/2] += ValV1[i];
00277     }
00278     EntropyV.Add(Entropy(ValV2));
00279     ValV1.MoveFrom(ValV2);
00280   }
00281   EntropyV.Reverse();
00282 }
00283 
00284 // solves for p: B = p * log2(p) + (1-p) * log2(1-p)
00285 double TSpecFunc::EntropyBias(const double& B){
00286   static TFltFltH BToP;
00287   if (BToP.Empty()) {
00288     for (double p = 0.5; p < 1.0; p +=0.0001) {
00289       double H = p * log(p) + (1.0-p) * log(1.0 - p);
00290       H = -H / log(2.0);
00291       BToP.AddDat(TMath::Round(H, 3), p);
00292     }
00293   }
00294   if (BToP.IsKey(TMath::Round(B, 3))) { return BToP.GetDat(TMath::Round(B, 3)); }
00295   else { return -1.0; }
00296 }
00297 
00298 // MLE of the power-law coefficient
00299 double TSpecFunc::GetPowerCoef(const TFltV& XValV, double MinX) {
00300   for (int i = 0; MinX <= 0.0 && i < XValV.Len(); i++) { 
00301     MinX = XValV[i]; }
00302   IAssert(MinX > 0.0);
00303   double LnSum=0.0;
00304   for (int i = 0; i < XValV.Len(); i++) {
00305     if (XValV[i].Val < MinX) continue;
00306     LnSum += log(XValV[i] / MinX);
00307   }
00308   return 1.0 + double(XValV.Len()) / LnSum;
00309 }
00310 
00311 double TSpecFunc::GetPowerCoef(const TFltPrV& XValCntV, double MinX) {
00312   for (int i = 0; MinX <= 0.0 && i < XValCntV.Len(); i++) { 
00313     MinX = XValCntV[i].Val1; }
00314   IAssert(MinX > 0.0);
00315   double NSamples=0.0, LnSum=0.0;
00316   for (int i = 0; i < XValCntV.Len(); i++) {
00317     if (XValCntV[i].Val1() < MinX) continue;
00318     LnSum += XValCntV[i].Val2 * log(XValCntV[i].Val1 / MinX);
00319     NSamples += XValCntV[i].Val2;
00320   }
00321   return 1.0 + NSamples / LnSum;
00322 }
00323 
00325 // Statistical-Moments
00326 TMom::TMom(const TFltV& _ValV):
00327   //WgtV(_ValV.Len(), 0), ValV(_ValV.Len(), 0),
00328   ValWgtV(_ValV.Len(), 0),
00329   SumW(), ValSumW(),
00330   UsableP(false), UnusableVal(-1),
00331   Mn(), Mx(),
00332   Mean(), Vari(), SDev(), SErr(),
00333   Median(), Quart1(), Quart3(),
00334   DecileV(), PercentileV(){
00335   for (int ValN=0; ValN<_ValV.Len(); ValN++){Add(_ValV[ValN], 1);}
00336   Def();
00337 }
00338 
00339 void TMom::Def(){
00340   IAssert(!DefP); DefP=true;
00341   UsableP=(SumW>0)&&(ValWgtV.Len()>0);
00342   if (UsableP){
00343     // Mn, Mx
00344     Mn=ValWgtV[0].Val1;
00345     Mx=ValWgtV[0].Val1;
00346     // Mean, Variance (Mn, Mx), Standard-Error
00347     Mean=ValSumW/SumW;
00348     Vari=0;
00349     if (ValWgtV.Len()>1){
00350       for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00351         const double Val=ValWgtV[ValN].Val1;
00352         Vari+=ValWgtV[ValN].Val2*TMath::Sqr(Val-Mean);
00353         if (Val<Mn){Mn=Val;}
00354         if (Val>Mx){Mx=Val;}
00355       }
00356       Vari=Vari/SumW;
00357       SErr=sqrt(Vari/(ValWgtV.Len()*(ValWgtV.Len()-1)));
00358     }
00359     // Standard-Deviation
00360     SDev=sqrt(double(Vari));
00361     // Median
00362     ValWgtV.Sort();
00363     double CurSumW = 0;
00364     for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00365       CurSumW += ValWgtV[ValN].Val2;
00366       if (CurSumW > 0.5*SumW) { 
00367         Median = ValWgtV[ValN].Val1; break; }
00368       else if (CurSumW == 0.5*SumW) {
00369         Median = 0.5 * (ValWgtV[ValN].Val1+ValWgtV[ValN+1].Val1); break; }
00370     }
00371     // Quartile-1 and Quartile-3
00372     Quart1=Quart3=TFlt::Mn;
00373     CurSumW = 0;
00374     for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00375       CurSumW += ValWgtV[ValN].Val2;
00376       if (Quart1==TFlt::Mn) {
00377         if (CurSumW > 0.25*SumW) {  Quart1 = ValWgtV[ValN].Val1; }
00378         //else if (CurSumW == 0.25*SumW) { Quart1 = 0.5 * (ValWgtV[ValN].Val1+ValWgtV[ValN+1].Val1); }
00379       } 
00380       if (Quart3==TFlt::Mn) {
00381         if (CurSumW > 0.75*SumW) { Quart3 = ValWgtV[ValN].Val1; }
00382         //else if (CurSumW == 0.75*SumW) { Quart3 = 0.5 * (ValWgtV[ValN].Val1+ValWgtV[ValN+1].Val1); }
00383       }
00384     }
00385     // Deciles & Percentiles
00386     CurSumW = 0;
00387     int DecileN = 1, PercentileN = 1;
00388     DecileV.Gen(11);  PercentileV.Gen(101);
00389     DecileV[0]=Mn; DecileV[10]=Mx;
00390     PercentileV[0]=Mn; PercentileV[100]=Mx;
00391     for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00392       CurSumW += ValWgtV[ValN].Val2;
00393       if (CurSumW > SumW*DecileN*0.1) { 
00394         DecileV[DecileN] = ValWgtV[ValN].Val1;  DecileN++; }
00395       if (CurSumW > SumW*PercentileN*0.01) {
00396         PercentileV[PercentileN] = ValWgtV[ValN].Val1;  PercentileN++; }
00397     }
00398   }
00399   ValWgtV.Clr();
00400 }
00401 
00402 
00403 
00404 double TMom::GetByNm(const TStr& MomNm) const {
00405   if (MomNm=="Mean"){return GetMean();}
00406   else if (MomNm=="Vari"){return GetVari();}
00407   else if (MomNm=="SDev"){return GetSDev();}
00408   else if (MomNm=="SErr"){return GetSErr();}
00409   else if (MomNm=="Median"){return GetMedian();}
00410   else if (MomNm=="Quart1"){return GetQuart1();}
00411   else if (MomNm=="Quart3"){return GetQuart3();}
00412   else if (MomNm=="Decile0"){return GetDecile(0);}
00413   else if (MomNm=="Decile1"){return GetDecile(1);}
00414   else if (MomNm=="Decile2"){return GetDecile(2);}
00415   else if (MomNm=="Decile3"){return GetDecile(3);}
00416   else if (MomNm=="Decile4"){return GetDecile(4);}
00417   else if (MomNm=="Decile5"){return GetDecile(5);}
00418   else if (MomNm=="Decile6"){return GetDecile(6);}
00419   else if (MomNm=="Decile7"){return GetDecile(7);}
00420   else if (MomNm=="Decile8"){return GetDecile(8);}
00421   else if (MomNm=="Decile9"){return GetDecile(9);}
00422   else if (MomNm=="Decile10"){return GetDecile(10);}
00423   else {Fail; return 0;}
00424 }
00425 
00426 TStr TMom::GetStrByNm(const TStr& MomNm, char* FmtStr) const {
00427   if (IsUsable()){
00428     if (FmtStr==NULL){
00429       return TFlt::GetStr(GetByNm(MomNm));
00430     } else {
00431       return TFlt::GetStr(GetByNm(MomNm), FmtStr);
00432     }
00433   } else {
00434     return "X";
00435   }
00436 }
00437 
00438 TStr TMom::GetStr(
00439  const char& SepCh, const char& DelimCh,
00440  const bool& DecileP, const bool& PercentileP, const TStr& FmtStr) const {
00441   TChA ChA;
00442   if (IsUsable()){
00443     ChA+="["; ChA+=SepCh;
00444     ChA+="Vals"; ChA+=DelimCh; ChA+=TInt::GetStr(GetVals()); ChA+=SepCh;
00445     ChA+="Min"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMn(), FmtStr.CStr()); ChA+=SepCh;
00446     ChA+="Max"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMx(), FmtStr.CStr()); ChA+=SepCh;
00447     ChA+="Mean"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMean(), FmtStr.CStr()); ChA+=SepCh;
00448     //ChA+="Vari"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetVari(), FmtStr.CStr()); ChA+=SepCh;
00449     ChA+="SDev"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetSDev(), FmtStr.CStr()); ChA+=SepCh;
00450     //ChA+="SErr"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetSErr(), FmtStr.CStr()); ChA+=SepCh;
00451     ChA+="Quart1"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetQuart1(), FmtStr.CStr()); ChA+=SepCh;
00452     ChA+="Median"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMedian(), FmtStr.CStr()); ChA+=SepCh;
00453     ChA+="Quart3"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetQuart3(), FmtStr.CStr()); ChA+=SepCh;
00454     if (DecileP){
00455       for (int DecileN=0; DecileN<=10; DecileN++){
00456         ChA+="Dec"; ChA+=TInt::GetStr(DecileN);
00457         ChA+=DelimCh; ChA+=TFlt::GetStr(GetDecile(DecileN), FmtStr.CStr());
00458         ChA+=SepCh;
00459       }
00460     }
00461     if (PercentileP){
00462       for (int PercentileN=0; PercentileN<=100; PercentileN++){
00463         ChA+="Per"; ChA+=TInt::GetStr(PercentileN);
00464         ChA+=DelimCh; ChA+=TFlt::GetStr(GetPercentile(PercentileN), FmtStr.CStr());
00465         ChA+=SepCh;
00466       }
00467     }
00468     ChA+="]";
00469   } else {
00470     ChA="[Unusable]";
00471   }
00472   return ChA;
00473 }
00474 
00475 TStr TMom::GetNmVStr(const TStr& VarPfx,
00476  const char& SepCh, const bool& DecileP, const bool& PercentileP){
00477   TChA ChA;
00478   ChA+=VarPfx; ChA+="Vals"; ChA+=SepCh;
00479   ChA+=VarPfx; ChA+="Min"; ChA+=SepCh;
00480   ChA+=VarPfx; ChA+="Max"; ChA+=SepCh;
00481   ChA+=VarPfx; ChA+="Mean"; ChA+=SepCh;
00482   //ChA+=VarPfx; ChA+="Vari"; ChA+=SepCh;
00483   ChA+=VarPfx; ChA+="SDev"; ChA+=SepCh;
00484   //ChA+=VarPfx; ChA+="SErr"; ChA+=SepCh;
00485   ChA+=VarPfx; ChA+="Quart1"; ChA+=SepCh;
00486   ChA+=VarPfx; ChA+="Median"; ChA+=SepCh;
00487   ChA+=VarPfx; ChA+="Quart3";
00488   if (DecileP){
00489     ChA+=SepCh;
00490     for (int DecileN=0; DecileN<=10; DecileN++){
00491       ChA+=VarPfx; ChA+="Dec"; ChA+=TInt::GetStr(DecileN);
00492       if (DecileN<10){ChA+=SepCh;}
00493     }
00494   }
00495   if (PercentileP){
00496     ChA+=SepCh;
00497     for (int PercentileN=0; PercentileN<=100; PercentileN++){
00498       ChA+=VarPfx; ChA+="Per"; ChA+=TInt::GetStr(PercentileN);
00499       if (PercentileN<100){ChA+=SepCh;}
00500     }
00501   }
00502   return ChA;
00503 }
00504 
00505 TStr TMom::GetValVStr(
00506  const char& SepCh, const bool& DecileP, const bool& PercentileP) const {
00507   TChA ChA;
00508   if (IsUsable()){
00509     ChA+=TInt::GetStr(GetVals()); ChA+=SepCh;
00510     ChA+=TFlt::GetStr(GetMn()); ChA+=SepCh;
00511     ChA+=TFlt::GetStr(GetMx()); ChA+=SepCh;
00512     ChA+=TFlt::GetStr(GetMean()); ChA+=SepCh;
00513     //ChA+=TFlt::GetStr(GetVari()); ChA+=SepCh;
00514     ChA+=TFlt::GetStr(GetSDev()); ChA+=SepCh;
00515     //ChA+=TFlt::GetStr(GetSErr()); ChA+=SepCh;
00516     ChA+=TFlt::GetStr(GetQuart1()); ChA+=SepCh;
00517     ChA+=TFlt::GetStr(GetMedian()); ChA+=SepCh;
00518     ChA+=TFlt::GetStr(GetQuart3()); ChA+=SepCh;
00519     if (DecileP){
00520       for (int DecileN=0; DecileN<=10; DecileN++){
00521         ChA+=TFlt::GetStr(GetDecile(DecileN)); ChA+=SepCh;
00522       }
00523     }
00524     if (PercentileP){
00525       for (int PercentileN=0; PercentileN<=100; PercentileN++){
00526         ChA+=TFlt::GetStr(GetPercentile(PercentileN)); ChA+=SepCh;
00527       }
00528     }
00529   } else {
00530     int Vals=8;
00531     if (DecileP){Vals+=11;}
00532     if (PercentileP){Vals+=101;}
00533     for (int ValN=0; ValN<Vals; ValN++){
00534       ChA="[Unusable]";
00535       if (ValN<Vals-1){ChA+=SepCh;}
00536     }
00537   }
00538   return ChA;
00539 }
00540 
00542 // Correlation
00543 TCorr::TCorr(const TFltV& ValV1, const TFltV& ValV2):
00544   ValVLen(ValV1.Len()), CorrCf(), CorrCfPrb(), FisherZ(){
00545   static const double TINY=1.0e-20;
00546   IAssert(ValV1.Len()==ValV2.Len());
00547 
00548   // calculate the means
00549   double MeanVal1=0; double MeanVal2=0;
00550   {for (int ValN=0; ValN<ValVLen; ValN++){
00551     MeanVal1+=ValV1[ValN]; MeanVal2+=ValV2[ValN];}}
00552   MeanVal1/=ValVLen; MeanVal2/=ValVLen;
00553 
00554   // calculate correlation coefficient
00555   double yt, xt;
00556   double syy=0.0; double sxy=0.0; double sxx=0.0;
00557   {for (int ValN=0; ValN<ValVLen; ValN++){
00558     xt=ValV1[ValN]-MeanVal1;
00559     yt=ValV2[ValN]-MeanVal2;
00560     sxx+=xt*xt;
00561     syy+=yt*yt;
00562     sxy+=xt*yt;
00563   }}
00564   if (sxx*syy==0){
00565     CorrCf=0; //** not in numerical recipes - check why (pojavi se, ko so same nicle)
00566   } else {
00567     CorrCf=sxy/sqrt(sxx*syy);
00568   }
00569   // calculate correlation coefficient significance level
00570   double df=ValVLen-2;
00571   double t=CorrCf*sqrt(df/((1.0-CorrCf+TINY)*(1.0+CorrCf+TINY)));
00572   CorrCfPrb=TSpecFunc::BetaI(0.5*df,0.5,df/(df+t*t));
00573   // calculate Fisher's Z transformation
00574   FisherZ=0.5*log((1.0+(CorrCf)+TINY)/(1.0-(CorrCf)+TINY));
00575 }
00576 
00578 // Statistical Tests
00579 void TStatTest::AveVar(const TFltV& ValV, double& Ave, double& Var){
00580   Ave=0;
00581   for (int ValN=0; ValN<ValV.Len(); ValN++){
00582     Ave+=ValV[ValN];}
00583   Ave/=ValV.Len();
00584   Var=0;
00585   double ep=0;
00586   for (int ValN=0; ValN<ValV.Len(); ValN++){
00587     double s=ValV[ValN]-Ave;
00588     ep+=s;
00589     Var+=s*s;
00590   }
00591   Var=(Var-ep*ep/ValV.Len())/(ValV.Len()-1);
00592 }
00593 
00594 double TStatTest::KsProb(const double& Alam) {
00595   const double EPS1 = 0.001;
00596   const double EPS2 = 1.0e-8;
00597   double a2 = -2.0*Alam*Alam, fac = 2.0, sum = 0.0, term, termbf = 0.0;
00598   for (int j=1; j <= 100; j++) {
00599     term = fac*exp(a2*j*j);
00600     sum += term;
00601     if (fabs(term) <= EPS1*termbf || fabs(term) <= EPS2*sum)
00602       return sum;
00603     fac = -fac;
00604     termbf = fabs(term);
00605   }
00606   return 1.0;
00607 }
00608 
00609 void TStatTest::ChiSquareOne(
00610  const TFltV& ObservedBinV, const TFltV& ExpectedBinV,
00611  double& ChiSquareVal, double& SignificancePrb){
00612   IAssert(ObservedBinV.Len()==ExpectedBinV.Len());
00613   int Bins=ObservedBinV.Len();
00614   int Constraints=0;
00615   int DegreesOfFreedom=Bins-Constraints;
00616   ChiSquareVal=0.0;
00617   for (int BinN=0; BinN<Bins; BinN++){
00618     IAssert(ExpectedBinV[BinN]>0);
00619     double BinDiff=ObservedBinV[BinN]-ExpectedBinV[BinN];
00620     ChiSquareVal+=BinDiff*BinDiff/ExpectedBinV[BinN];
00621   }
00622   SignificancePrb=
00623    TSpecFunc::GammaQ(0.5*(DegreesOfFreedom), 0.5*(ChiSquareVal));
00624 }
00625 
00626 void TStatTest::ChiSquareTwo(
00627  const TFltV& ObservedBin1V, const TFltV& ObservedBin2V,
00628  double& ChiSquareVal, double& SignificancePrb){
00629   IAssert(ObservedBin1V.Len()==ObservedBin1V.Len());
00630   int Bins=ObservedBin1V.Len();
00631   int Constraints=0;
00632   int DegreesOfFreedom=Bins-Constraints;
00633   ChiSquareVal=0.0;
00634   for (int BinN=0; BinN<Bins; BinN++){
00635     if ((ObservedBin1V[BinN]==0.0) && (ObservedBin2V[BinN]==0.0)){
00636       DegreesOfFreedom--;
00637     } else {
00638       double BinDiff=ObservedBin1V[BinN]-ObservedBin2V[BinN];
00639       ChiSquareVal+=BinDiff*BinDiff/(ObservedBin1V[BinN]+ObservedBin2V[BinN]);
00640     }
00641   }
00642   SignificancePrb=
00643    TSpecFunc::GammaQ(0.5*(DegreesOfFreedom),0.5*(ChiSquareVal));
00644 }
00645 
00646 void TStatTest::TTest(
00647  const TFltV& ValV1, const TFltV& ValV2, double& TTestVal, double& TTestPrb){
00648   /*double Ave1; double Var1;
00649   AveVar(ValV1, Ave1, Var1);
00650   double Ave2; double Var2;
00651   AveVar(ValV2, Ave2, Var2);*/
00652 
00653   PMom Val1Mom=TMom::New(ValV1);
00654   PMom Val2Mom=TMom::New(ValV2);
00655   double ave1=Val1Mom->GetMean();
00656   double ave2=Val2Mom->GetMean();
00657   double var1=Val1Mom->GetVari();
00658   double var2=Val2Mom->GetVari();
00659   int n1=ValV1.Len();
00660   int n2=ValV2.Len();
00661 
00662   TTestVal=(ave1-ave2)/sqrt(var1/n1+var2/n2);
00663   double df=TMath::Sqr(var1/n1+var2/n2)/(TMath::Sqr(var1/n1)/(n1-1)+TMath::Sqr(var2/n2)/(n2-1));
00664   TTestPrb=TSpecFunc::BetaI(0.5*df, 0.5, df/(df+TMath::Sqr(TTestVal)));
00665 }
00666 
00667 void TStatTest::KsTest(const TFltV& ValV1, const TFltV& ValV2, double& DStat, double& PVal) {
00668   IAssert(ValV1.IsSorted() && ValV2.IsSorted());
00669   int i1=0, i2=0;
00670   double CumSum1=0.0, CumSum2=0.0, Cdf1=0.0, Cdf2=0.0;
00671   const double N1 = ValV1.Len();
00672   const double N2 = ValV2.Len();
00673   if (! (N1 > 0.0 && N2 > 0.0)) { DStat = 1.0;  PVal = 0.0;  return; }
00674   DStat=0.0; PVal=0.0;
00675   while (i1 < ValV1.Len() && i2 < ValV2.Len()) {
00676     const double X1 = ValV1[i1];
00677     const double X2 = ValV2[i2];
00678     if (X1 <= X2) {
00679       CumSum1 += 1;
00680       Cdf1 = (CumSum1 / N1);
00681       i1++;
00682     }
00683     if (X2 <= X1) {
00684       CumSum2 += 1;
00685       Cdf2 = (CumSum2 / N2);
00686       i2++;
00687     }
00688     DStat = TMath::Mx(DStat, fabs(Cdf1 - Cdf2));
00689   }
00690   const double En = sqrt( N1*N2 / (N1+N2));
00691   PVal = TStatTest::KsProb((En+0.12+0.11/En)*DStat);
00692 }
00693 
00694 void TStatTest::KsTest(const TFltPrV& ValCntV1, const TFltPrV& ValCntV2, double& DStat, double& PVal) {
00695   IAssert(ValCntV1.IsSorted() && ValCntV2.IsSorted());
00696   int i1=0, i2=0;
00697   double N1=0.0, N2=0.0, CumSum1=0.0, CumSum2=0.0, Cdf1=0.0, Cdf2=0.0;
00698   DStat=0.0;  PVal=0.0;
00699   for (int i = 0; i < ValCntV1.Len(); i++) N1 += ValCntV1[i].Val2;
00700   for (int i = 0; i < ValCntV2.Len(); i++) N2 += ValCntV2[i].Val2;
00701   if (! (N1 > 0.0 && N2 > 0.0)) { DStat = 1.0;  PVal = 0.0;  return; }
00702 
00703   while (i1 < ValCntV1.Len() && i2 < ValCntV2.Len()) {
00704     const double X1 = ValCntV1[i1].Val1;
00705     const double X2 = ValCntV2[i2].Val1;
00706     if (X1 <= X2) {
00707       CumSum1 += ValCntV1[i1].Val2;
00708       Cdf1 = (CumSum1 / N1);
00709       i1++;
00710     }
00711     if (X2 <= X1) {
00712       CumSum2 += ValCntV2[i2].Val2;
00713       Cdf2 = (CumSum2 / N2);
00714       i2++;
00715     }
00716     DStat = TMath::Mx(DStat, fabs(Cdf1 - Cdf2));
00717   }
00718   const double En = sqrt( N1*N2 / (N1+N2));
00719   PVal = TStatTest::KsProb((En+0.12+0.11/En)*DStat);
00720 }
00721 
00723 // Combinations
00724 bool TComb::GetNext(){
00725   if (ItemV.Len()==0){
00726     ItemV.Gen(Order, Order);
00727     for (int OrderN=0; OrderN<Order; OrderN++){
00728       ItemV[OrderN]=OrderN;}
00729     return true;
00730   } else {
00731     if (ItemV.Last()==Items-1){
00732       int OrderN=Order-1;
00733       while ((OrderN>=0)&&(ItemV[OrderN]==Items-(Order-OrderN-1)-1)){OrderN--;}
00734       if (OrderN<0){
00735         return false;
00736       } else {
00737         ItemV[OrderN]++;
00738         for (int SubOrderN=OrderN+1; SubOrderN<Order; SubOrderN++){
00739           ItemV[SubOrderN]=ItemV[SubOrderN-1]+1;}
00740         CombN++; return true;
00741       }
00742     } else {
00743       ItemV.Last()++; CombN++; return true;
00744     }
00745   }
00746 }
00747 
00748 int TComb::GetCombs() const {
00749   int LCombs=1; int HCombs=1;
00750   for (int OrderN=0; OrderN<Order; OrderN++){
00751     LCombs*=OrderN+1; HCombs*=Items-OrderN;}
00752   int Combs=HCombs/LCombs;
00753   return Combs;
00754 }
00755 
00756 void TComb::Wr(){
00757   printf("%d:[", GetCombN());
00758   for (int OrderN=0; OrderN<Order; OrderN++){
00759     if (OrderN>0){printf(" ");}
00760     printf("%d", ItemV[OrderN]());
00761   }
00762   printf("]\n");
00763 }
00764 
00766 // Linear-Regression
00767 PLinReg TLinReg::New(const TFltVV& _XVV, const TFltV& _YV, const TFltV& _SigV){
00768   PLinReg LinReg=PLinReg(new TLinReg());
00769   LinReg->XVV=_XVV;
00770   LinReg->YV=_YV;
00771   if (_SigV.Empty()){
00772     LinReg->SigV.Gen(LinReg->YV.Len());
00773     LinReg->SigV.PutAll(1);
00774   } else {
00775     LinReg->SigV=_SigV;
00776   }
00777   LinReg->Recs=LinReg->XVV.GetXDim();
00778   LinReg->Vars=LinReg->XVV.GetYDim();
00779   IAssert(LinReg->Recs>0);
00780   IAssert(LinReg->Vars>0);
00781   IAssert(LinReg->YV.Len()==LinReg->Recs);
00782   IAssert(LinReg->SigV.Len()==LinReg->Recs);
00783   LinReg->CovarVV.Gen(LinReg->Vars+1, LinReg->Vars+1);
00784   LinReg->CfV.Gen(LinReg->Vars+1);
00785   LinReg->NR_lfit();
00786   return LinReg;
00787 }
00788 
00789 void TLinReg::NR_covsrt(
00790  TFltVV& CovarVV, const int& Vars, const TIntV& ia, const int& mfit){
00791   for (int i=mfit+1; i<=Vars; i++){
00792     for (int j=1; j<=i; j++){
00793       CovarVV.At(i, j)=0; CovarVV.At(j, i)=0.0;}
00794   }
00795   int k=mfit;
00796   for (int j=Vars; j>=1; j--){
00797     if (ia[j]!=0){
00798       for (int i=1; i<=Vars; i++){Swap(CovarVV.At(i, k), CovarVV.At(i, j));}
00799       {for (int i=1; i<=Vars; i++){Swap(CovarVV.At(k, i), CovarVV.At(j, i));}}
00800       k--;
00801     }
00802   }
00803 }
00804 
00805 void TLinReg::NR_gaussj(TFltVV& a, const int& n, TFltVV& b, const int& m){
00806   int i, icol, irow=0, j, k, l, ll;
00807   double big, dum, pivinv;
00808 
00809   TIntV indxc(n+1);
00810   TIntV indxr(n+1);
00811   TIntV ipiv(n+1);
00812   for (j=1; j<=n; j++){ipiv[j]=0;}
00813   for (i=1; i<=n; i++){
00814     big=0.0;
00815     for (j=1; j<=n; j++){
00816       if (ipiv[j]!=1){
00817         for (k=1; k<=n; k++){
00818           if (ipiv[k]==0){
00819             if (fabs(double(a.At(j, k))) >= big){
00820               big=fabs(double(a.At(j, k)));
00821               irow=j;
00822               icol=k;
00823             }
00824           } else
00825           if (ipiv[k]>1){
00826             TExcept::Throw("Singular Matrix(1) in Gauss");}
00827         }
00828       }
00829     }
00830     ipiv[icol]++;
00831     if (irow != icol){
00832       for (l=1; l<=n; l++){Swap(a.At(irow, l), a.At(icol, l));}
00833       for (l=1; l<=m; l++){Swap(b.At(irow, l), b.At(icol, l));}
00834     }
00835     indxr[i]=irow;
00836     indxc[i]=icol;
00837     if (a.At(icol, icol)==0.0){
00838       TExcept::Throw("Singular Matrix(1) in Gauss");}
00839     pivinv=1.0/a.At(icol, icol);
00840     a.At(icol, icol)=1.0;
00841     for (l=1; l<=n; l++){a.At(icol, l)=a.At(icol, l)*pivinv;}
00842     for (l=1; l<=m; l++){b.At(icol, l)=b.At(icol, l)*pivinv;}
00843     for (ll=1; ll<=n; ll++){
00844       if (ll != icol){
00845         dum=a.At(ll, icol);
00846         a.At(ll, icol)=0.0;
00847         for (l=1;l<=n;l++){a.At(ll, l)-=a.At(icol, l)*dum;}
00848         for (l=1;l<=m;l++){b.At(ll, l)-=b.At(icol, l)*dum;}
00849       }
00850     }
00851   }
00852   for (l=n; l>=1; l--){
00853     if (indxr[l]!=indxc[l]){
00854       for (k=1; k<=n; k++){
00855         Swap(a.At(k, indxr[l]), a.At(k, indxc[l]));}
00856     }
00857   }
00858 }
00859 
00860 void TLinReg::NR_lfit(){
00861   int i,j,k,l,m,mfit=0;
00862   double ym,wt,sum,sig2i;
00863 
00864   TIntV ia(Vars+1); for (i=1; i<=Vars; i++){ia[i]=1;}
00865   TFltVV beta(Vars+1, 1+1);
00866   TFltV afunc(Vars+1);
00867   for (j=1;j<=Vars;j++){
00868     if (ia[j]!=0){mfit++;}}
00869   if (mfit==0){TExcept::Throw("No parameters to be fitted in LFit");}
00870   for (j=1; j<=mfit; j++){
00871     for (k=1; k<=mfit; k++){CovarVV.At(j, k)=0.0;}
00872     beta.At(j, 1)=0.0;
00873   }
00874   for (i=1; i<=Recs; i++){
00875     GetXV(i, afunc); // funcs(XVV[i],afunc,Vars);
00876     ym=GetY(i);
00877     if (mfit<Vars){
00878       for (j=1;j<=Vars;j++){
00879         if (ia[j]==0){ym-=CfV[j]*afunc[j];}}
00880     }
00881     sig2i=1.0/TMath::Sqr(GetSig(i));
00882     for (j=0, l=1; l<=Vars; l++){
00883       if (ia[l]!=0){
00884         wt=afunc[l]*sig2i;
00885         for (j++, k=0, m=1; m<=l; m++){
00886           if (ia[m]!=0){CovarVV.At(j, ++k)+=wt*afunc[m];}
00887         }
00888         beta.At(j, 1)+=ym*wt;
00889       }
00890     }
00891   }
00892   for (j=2; j<=mfit; j++){
00893     for (k=1; k<j; k++){CovarVV.At(k, j)=CovarVV.At(j, k);}
00894   }
00895   NR_gaussj(CovarVV, mfit, beta, 1);
00896   for (j=0, l=1; l<=Vars; l++){
00897     if (ia[l]!=0){CfV[l]=beta.At(++j, 1);}
00898   }
00899   ChiSq=0.0;
00900   for (i=1; i<=Recs; i++){
00901     GetXV(i, afunc); // funcs(XVV[i],afunc,Vars);
00902     for (sum=0.0, j=1; j<=Vars; j++){sum+=CfV[j]*afunc[j];}
00903     ChiSq+=TMath::Sqr((GetY(i)-sum)/GetSig(i));
00904   }
00905   NR_covsrt(CovarVV, Vars, ia, mfit);
00906 }
00907 
00908 void TLinReg::Wr() const {
00909   printf("\n%11s %21s\n","parameter","uncertainty");
00910   for (int i=0; i<Vars;i++){
00911     printf("  a[%1d] = %8.6f %12.6f\n",
00912      i+1, GetCf(i), GetCfUncer(i));
00913   }
00914   printf("chi-squared = %12f\n", GetChiSq());
00915 
00916   printf("full covariance matrix\n");
00917   {for (int i=0;i<Vars;i++){
00918     for (int j=0;j<Vars;j++){
00919       printf("%12f", GetCovar(i, j));}
00920     printf("\n");
00921   }}
00922 }
00923 
00925 // Singular-Value-Decomposition
00926 PSvd TSvd::New(const TFltVV& _XVV, const TFltV& _YV, const TFltV& _SigV){
00927   PSvd Svd=PSvd(new TSvd());
00928   Svd->XVV=_XVV;
00929   Svd->YV=_YV;
00930   if (_SigV.Empty()){
00931     Svd->SigV.Gen(Svd->YV.Len());
00932     Svd->SigV.PutAll(1);
00933   } else {
00934     Svd->SigV=_SigV;
00935   }
00936   Svd->Recs=Svd->XVV.GetXDim();
00937   Svd->Vars=Svd->XVV.GetYDim();
00938   IAssert(Svd->Recs>0);
00939   IAssert(Svd->Vars>0);
00940   IAssert(Svd->YV.Len()==Svd->Recs);
00941   IAssert(Svd->SigV.Len()==Svd->Recs);
00942   Svd->CovarVV.Gen(Svd->Vars+1, Svd->Vars+1);
00943   Svd->CfV.Gen(Svd->Vars+1);
00944   Svd->NR_svdfit();
00945   return Svd;
00946 }
00947 
00948 double TSvd::NR_pythag(double a, double b){
00949   double absa,absb;
00950   absa=fabs(a);
00951   absb=fabs(b);
00952   if (absa > absb){
00953     return absa*sqrt(1.0+TMath::Sqr(absb/absa));
00954   } else {
00955     return (absb == 0.0 ? 0.0 : absb*sqrt(1.0+TMath::Sqr(absa/absb)));
00956   }
00957 }
00958 
00959 void TSvd::NR_svdcmp(TFltVV& a, int m, int n, TFltV& w, TFltVV& v){
00960   int flag,i,its,j,jj,k,l=0,nm;
00961   double anorm,c,f,g,h,s,scale,x,y,z;
00962 
00963   TFltV rv1(n+1);
00964   g=scale=anorm=0.0;
00965   for (i=1;i<=n;i++) {
00966     l=i+1;
00967     rv1[i]=scale*g;
00968     g=s=scale=0.0;
00969     if (i <= m) {
00970       for (k=i;k<=m;k++) scale += fabs(double(a.At(k,i)));
00971       if (scale) {
00972         for (k=i;k<=m;k++) {
00973           a.At(k,i) /= scale;
00974           s += a.At(k,i)*a.At(k,i);
00975         }
00976         f=a.At(i,i);
00977         g = -NR_SIGN(sqrt(s),f);
00978         h=f*g-s;
00979         a.At(i,i)=f-g;
00980         for (j=l;j<=n;j++) {
00981           for (s=0.0,k=i;k<=m;k++) s += a.At(k,i)*a(k,j);
00982           f=s/h;
00983           for (k=i;k<=m;k++) a.At(k,j) += f*a.At(k,i);
00984         }
00985         for (k=i;k<=m;k++) a.At(k,i) *= scale;
00986       }
00987     }
00988     w[i]=scale *g;
00989     g=s=scale=0.0;
00990     if (i <= m && i != n) {
00991       for (k=l;k<=n;k++) scale += fabs(double(a.At(i,k)));
00992       if (scale) {
00993         for (k=l;k<=n;k++) {
00994           a.At(i,k) /= scale;
00995           s += a.At(i,k)*a.At(i,k);
00996         }
00997         f=a.At(i,l);
00998         g = -NR_SIGN(sqrt(s),f);
00999         h=f*g-s;
01000         a.At(i,l)=f-g;
01001         for (k=l;k<=n;k++) rv1[k]=a.At(i,k)/h;
01002         for (j=l;j<=m;j++) {
01003           for (s=0.0,k=l;k<=n;k++) s += a.At(j,k)*a.At(i,k);
01004           for (k=l;k<=n;k++) a.At(j,k) += s*rv1[k];
01005         }
01006         for (k=l;k<=n;k++) a.At(i,k) *= scale;
01007       }
01008     }
01009     anorm=NR_FMAX(anorm,(fabs(double(w[i]))+fabs(double(rv1[i]))));
01010   }
01011   for (i=n;i>=1;i--) {
01012     if (i < n) {
01013       if (g) {
01014         for (j=l;j<=n;j++)
01015           v.At(j,i)=(a.At(i,j)/a.At(i,l))/g;
01016         for (j=l;j<=n;j++) {
01017           for (s=0.0,k=l;k<=n;k++) s += a.At(i,k)*v.At(k,j);
01018           for (k=l;k<=n;k++) v.At(k,j) += s*v.At(k,i);
01019         }
01020       }
01021       for (j=l;j<=n;j++) v.At(i,j)=v.At(j,i)=0.0;
01022     }
01023     v.At(i,i)=1.0;
01024     g=rv1[i];
01025     l=i;
01026   }
01027   for (i=NR_IMIN(m,n);i>=1;i--) {
01028     l=i+1;
01029     g=w[i];
01030     for (j=l;j<=n;j++) a.At(i,j)=0.0;
01031     if (g) {
01032       g=1.0/g;
01033       for (j=l;j<=n;j++) {
01034         for (s=0.0,k=l;k<=m;k++) s += a.At(k,i)*a.At(k,j);
01035         f=(s/a.At(i,i))*g;
01036         for (k=i;k<=m;k++) a.At(k,j) += f*a.At(k,i);
01037       }
01038       for (j=i;j<=m;j++) a.At(j,i) *= g;
01039     } else for (j=i;j<=m;j++) a.At(j,i)=0.0;
01040     a.At(i,i)++;
01041   }
01042   for (k=n;k>=1;k--) {
01043     for (its=1;its<=30;its++) {
01044       flag=1;
01045       for (l=k;l>=1;l--) {
01046         nm=l-1;
01047         if ((double)(fabs(double(rv1[l])+anorm)) == anorm) {
01048           flag=0;
01049           break;
01050         }
01051         if ((double)(fabs(double(w[nm]))+anorm) == anorm) break;
01052       }
01053       if (flag) {
01054         c=0.0;
01055         s=1.0;
01056         for (i=l;i<=k;i++) {
01057           f=s*rv1[i];
01058           rv1[i]=c*rv1[i];
01059           if ((double)(fabs(f)+anorm) == anorm) break;
01060           g=w[i];
01061           h=NR_pythag(f,g);
01062           w[i]=h;
01063           h=1.0/h;
01064           c=g*h;
01065           s = -f*h;
01066           for (j=1;j<=m;j++) {
01067             y=a.At(j,nm);
01068             z=a.At(j,i);
01069             a.At(j,nm)=y*c+z*s;
01070             a.At(j,i)=z*c-y*s;
01071           }
01072         }
01073       }
01074       z=w[k];
01075       if (l == k) {
01076         if (z < 0.0) {
01077           w[k] = -z;
01078           for (j=1;j<=n;j++) v.At(j,k) = -v.At(j,k);
01079         }
01080         break;
01081       }
01082       if (its==30){
01083         TExcept::Throw("no convergence in 30 svdcmp iterations");}
01084       x=w[l];
01085       nm=k-1;
01086       y=w[nm];
01087       g=rv1[nm];
01088       h=rv1[k];
01089       f=((y-z)*(y+z)+(g-h)*(g+h))/(2.0*h*y);
01090       g=NR_pythag(f,1.0);
01091       f=((x-z)*(x+z)+h*((y/(f+NR_SIGN(g,f)))-h))/x;
01092       c=s=1.0;
01093       for (j=l;j<=nm;j++) {
01094         i=j+1;
01095         g=rv1[i];
01096         y=w[i];
01097         h=s*g;
01098         g=c*g;
01099         z=NR_pythag(f,h);
01100         rv1[j]=z;
01101         c=f/z;
01102         s=h/z;
01103         f=x*c+g*s;
01104         g = g*c-x*s;
01105         h=y*s;
01106         y *= c;
01107         for (jj=1;jj<=n;jj++) {
01108           x=v.At(jj,j);
01109           z=v.At(jj,i);
01110           v.At(jj,j)=x*c+z*s;
01111           v.At(jj,i)=z*c-x*s;
01112         }
01113         z=NR_pythag(f,h);
01114         w[j]=z;
01115         if (z) {
01116           z=1.0/z;
01117           c=f*z;
01118           s=h*z;
01119         }
01120         f=c*g+s*y;
01121         x=c*y-s*g;
01122         for (jj=1;jj<=m;jj++) {
01123           y=a.At(jj,j);
01124           z=a.At(jj,i);
01125           a.At(jj,j)=y*c+z*s;
01126           a.At(jj,i)=z*c-y*s;
01127         }
01128       }
01129       rv1[l]=0.0;
01130       rv1[k]=f;
01131       w[k]=x;
01132     }
01133   }
01134 }
01135 
01136 void TSvd::NR_svbksb(
01137  TFltVV& u, TFltV& w, TFltVV& v, int m, int n, TFltV& b, TFltV& x){
01138   int jj,j,i;
01139   double s;
01140 
01141   TFltV tmp(n+1);
01142   for (j=1;j<=n;j++) {
01143     s=0.0;
01144     if (w[j]) {
01145       for (i=1;i<=m;i++) s += u.At(i,j)*b[i];
01146       s /= w[j];
01147     }
01148     tmp[j]=s;
01149   }
01150   for (j=1;j<=n;j++) {
01151     s=0.0;
01152     for (jj=1;jj<=n;jj++) s += v.At(j,jj)*tmp[jj];
01153     x[j]=s;
01154   }
01155 }
01156 
01157 void TSvd::NR_svdvar(TFltVV& v, int ma, TFltV& w, TFltVV& cvm){
01158   int k,j,i;
01159   double sum;
01160 
01161   TFltV wti(ma+1);
01162   for (i=1;i<=ma;i++) {
01163     wti[i]=0.0;
01164     if (w[i]) wti[i]=1.0/(w[i]*w[i]);
01165   }
01166   for (i=1;i<=ma;i++) {
01167     for (j=1;j<=i;j++) {
01168       for (sum=0.0,k=1;k<=ma;k++) sum += v.At(i,k)*v.At(j,k)*wti[k];
01169       cvm.At(j,i)=cvm.At(i,j)=sum;
01170     }
01171   }
01172 }
01173 
01174 void TSvd::NR_svdfit(){
01175   int j,i;
01176   double wmax,tmp,thresh,sum;
01177   double TOL=1.0e-5;
01178 
01179   TFltVV u(Recs+1, Vars+1);
01180   TFltVV v(Vars+1, Vars+1);
01181   TFltV w(Vars+1);
01182   TFltV b(Recs+1);
01183   TFltV afunc(Vars+1);
01184   for (i=1;i<=Recs;i++) {
01185     GetXV(i, afunc); // (*funcs)(x[i],afunc,Vars);
01186     tmp=1.0/GetSig(i);
01187     for (j=1;j<=Vars;j++){u.At(i,j)=afunc[j]*tmp;}
01188     b[i]=GetY(i)*tmp;
01189   }
01190   NR_svdcmp(u,Recs,Vars,w,v);
01191   wmax=0.0;
01192   for (j=1;j<=Vars;j++){
01193     if (w[j] > wmax){wmax=w[j];}}
01194   thresh=TOL*wmax;
01195   for (j=1;j<=Vars;j++){
01196     if (double(w[j])<thresh){w[j]=0.0;}}
01197   NR_svbksb(u,w,v,Recs,Vars,b,CfV);
01198   ChiSq=0.0;
01199   for (i=1;i<=Recs;i++) {
01200     GetXV(i, afunc); // (*funcs)(x[i],afunc,Vars);
01201     for (sum=0.0,j=1;j<=Vars;j++){sum += CfV[j]*afunc[j];}
01202     ChiSq += (tmp=(GetY(i)-sum)/GetSig(i),tmp*tmp);
01203   }
01204   // covariance matrix calculation
01205   CovarVV.Gen(Vars+1, Vars+1);
01206   NR_svdvar(v, Vars, w, CovarVV);
01207 }
01208 
01209 void TSvd::GetCfV(TFltV& _CfV){
01210   _CfV=CfV; _CfV.Del(0);
01211 }
01212 
01213 void TSvd::GetCfUncerV(TFltV& CfUncerV){
01214   CfUncerV.Gen(Vars);
01215   for (int VarN=0; VarN<Vars; VarN++){
01216     CfUncerV[VarN]=GetCfUncer(VarN);
01217   }
01218 }
01219 
01220 // all vectors (matrices) start with 0 
01221 void TSvd::Svd(const TFltVV& InMtx, TFltVV& LSingV, TFltV& SingValV, TFltVV& RSingV) {
01222   //LSingV = InMtx;
01223   LSingV.Gen(InMtx.GetYDim()+1, InMtx.GetYDim()+1);
01224   // create 1 based adjacency matrix
01225   for (int x = 0; x < InMtx.GetXDim(); x++) {
01226     for (int y = 0; y < InMtx.GetYDim(); y++) {
01227       LSingV.At(x+1, y+1) = InMtx.At(x, y);
01228     }
01229   }
01230   RSingV.Gen(InMtx.GetYDim()+1, InMtx.GetYDim()+1);
01231   SingValV.Gen(InMtx.GetYDim()+1);
01232   TSvd::NR_svdcmp(LSingV, InMtx.GetXDim(), InMtx.GetYDim(), SingValV, RSingV);
01233   // 0-th singular value/vector is full of zeros, delete it
01234   SingValV.Del(0);
01235   LSingV.DelX(0); LSingV.DelY(0);
01236   RSingV.DelX(0); RSingV.DelY(0);
01237 }
01238 
01239 // InMtx1 is 1-based (0-th row/column are full of zeros)
01240 // returned singular vectors/values are 0 based
01241 void TSvd::Svd1Based(const TFltVV& InMtx1, TFltVV& LSingV, TFltV& SingValV, TFltVV& RSingV) {
01242   LSingV = InMtx1;
01243   SingValV.Gen(InMtx1.GetYDim());
01244   RSingV.Gen(InMtx1.GetYDim(), InMtx1.GetYDim());
01245   TSvd::NR_svdcmp(LSingV, InMtx1.GetXDim()-1, InMtx1.GetYDim()-1, SingValV, RSingV);
01246   // 0-th singular value/vector is full of zeros, delete it
01247   SingValV.Del(0);
01248   LSingV.DelX(0); LSingV.DelY(0);
01249   RSingV.DelX(0); RSingV.DelY(0);
01250 }
01251 
01252 void TSvd::Wr() const {
01253   printf("\n%11s %21s\n","parameter","uncertainty");
01254   for (int i=0; i<Vars;i++){
01255     printf("  a[%1d] = %8.6f %12.6f\n",
01256      i+1, GetCf(i), GetCfUncer(i));
01257   }
01258   printf("chi-squared = %12f\n", GetChiSq());
01259 
01260   printf("full covariance matrix\n");
01261   {for (int i=0;i<Vars;i++){
01262     for (int j=0;j<Vars;j++){
01263       printf("%12f", GetCovar(i, j));}
01264     printf("\n");
01265   }}
01266 }
01267 
01269 // Histogram
01270 void THist::Add(const double& Val, const bool& OnlyInP) {
01271         // get bucket number
01272     const int BucketN = int(floor((Val - MnVal) / BucketSize));
01273         if (OnlyInP) { 
01274                 // value should be inside
01275                 EAssert(MnVal <= Val && Val <= MxVal);
01276                 BucketV[BucketN]++;
01277         } else {
01278                 // value does not need to be inside
01279                 if (BucketN < 0) {
01280                         BucketV[0]++;
01281                 } else if (BucketN < BucketV.Len()) {
01282                         BucketV[BucketN]++;
01283                 } else {
01284                         BucketV.Last()++;
01285                 }
01286         }
01287         // for computing percentage
01288         Vals++;
01289 }
01290 
01291 void THist::SaveStat(const TStr& ValNm, TSOut& FOut) const {
01292     FOut.PutStrLn("#" + ValNm + ": " + Vals.GetStr());
01293     const int Buckets = BucketV.Len() - 1;
01294     for (int BucketN = 0; BucketN < Buckets; BucketN++) {
01295         FOut.PutStrLn(TStr::Fmt("%d-%d\t%d", BucketSize*BucketN,
01296             BucketSize*(BucketN+1), BucketV[BucketN]()));
01297     }
01298     if (BucketV.Last() > 0) {
01299         FOut.PutStrLn(TStr::Fmt("%d-\t%d", BucketSize*Buckets, BucketV.Last()()));
01300     }
01301 
01302 }