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
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))); //J: This seems to be wrong
00358       if (Vari > 0.0 && SumW > 0.0) {
00359         SErr=sqrt(double(Vari))/sqrt(double(SumW)); //J: This seems to ok: StdDev/sqrt(N)
00360       } else { SErr = Mx; } // some big number
00361     }
00362     // Standard-Deviation
00363     SDev=sqrt(double(Vari));
00364     // Median
00365     ValWgtV.Sort();
00366     double CurSumW = 0;
00367     for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00368       CurSumW += ValWgtV[ValN].Val2;
00369       if (CurSumW > 0.5*SumW) { 
00370         Median = ValWgtV[ValN].Val1; break; }
00371       else if (CurSumW == 0.5*SumW) {
00372         Median = 0.5 * (ValWgtV[ValN].Val1+ValWgtV[ValN+1].Val1); break; }
00373     }
00374     // Quartile-1 and Quartile-3
00375     Quart1=Quart3=TFlt::Mn;
00376     CurSumW = 0;
00377     for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00378       CurSumW += ValWgtV[ValN].Val2;
00379       if (Quart1==TFlt::Mn) {
00380         if (CurSumW > 0.25*SumW) {  Quart1 = ValWgtV[ValN].Val1; }
00381         //else if (CurSumW == 0.25*SumW) { Quart1 = 0.5 * (ValWgtV[ValN].Val1+ValWgtV[ValN+1].Val1); }
00382       } 
00383       if (Quart3==TFlt::Mn) {
00384         if (CurSumW > 0.75*SumW) { Quart3 = ValWgtV[ValN].Val1; }
00385         //else if (CurSumW == 0.75*SumW) { Quart3 = 0.5 * (ValWgtV[ValN].Val1+ValWgtV[ValN+1].Val1); }
00386       }
00387     }
00388     // Mode (value with max total weight)
00389     THash<TFlt, TFlt> ValWgtH;
00390     for (int i = 0; i < ValWgtV.Len(); i++) {
00391       ValWgtH.AddDat(ValWgtV[i].Val1) += ValWgtV[i].Val2; }
00392     Mode = TFlt::Mn; double MxWgt=TFlt::Mn;
00393     for (int v = 0; v < ValWgtH.Len(); v++) {
00394       if (ValWgtH[v] > MxWgt) { MxWgt=ValWgtH[v]; Mode=ValWgtH.GetKey(v); }
00395     }
00396     // Deciles & Percentiles
00397     CurSumW = 0;
00398     int DecileN = 1, PercentileN = 1;
00399     DecileV.Gen(11);  PercentileV.Gen(101);
00400     DecileV[0]=Mn; DecileV[10]=Mx;
00401     PercentileV[0]=Mn; PercentileV[100]=Mx;
00402     for (int ValN=0; ValN<ValWgtV.Len(); ValN++){
00403       CurSumW += ValWgtV[ValN].Val2;
00404       if (CurSumW > SumW*DecileN*0.1) { 
00405         DecileV[DecileN] = ValWgtV[ValN].Val1;  DecileN++; }
00406       if (CurSumW > SumW*PercentileN*0.01) {
00407         PercentileV[PercentileN] = ValWgtV[ValN].Val1;  PercentileN++; }
00408     }
00409   }
00410   ValWgtV.Clr();
00411 }
00412 
00413 
00414 
00415 double TMom::GetByNm(const TStr& MomNm) const {
00416   if (MomNm=="Mean"){return GetMean();}
00417   else if (MomNm=="Vari"){return GetVari();}
00418   else if (MomNm=="SDev"){return GetSDev();}
00419   else if (MomNm=="SErr"){return GetSErr();}
00420   else if (MomNm=="Median"){return GetMedian();}
00421   else if (MomNm=="Quart1"){return GetQuart1();}
00422   else if (MomNm=="Quart3"){return GetQuart3();}
00423   else if (MomNm=="Decile0"){return GetDecile(0);}
00424   else if (MomNm=="Decile1"){return GetDecile(1);}
00425   else if (MomNm=="Decile2"){return GetDecile(2);}
00426   else if (MomNm=="Decile3"){return GetDecile(3);}
00427   else if (MomNm=="Decile4"){return GetDecile(4);}
00428   else if (MomNm=="Decile5"){return GetDecile(5);}
00429   else if (MomNm=="Decile6"){return GetDecile(6);}
00430   else if (MomNm=="Decile7"){return GetDecile(7);}
00431   else if (MomNm=="Decile8"){return GetDecile(8);}
00432   else if (MomNm=="Decile9"){return GetDecile(9);}
00433   else if (MomNm=="Decile10"){return GetDecile(10);}
00434   else {Fail; return 0;}
00435 }
00436 
00437 TStr TMom::GetStrByNm(const TStr& MomNm, char* FmtStr) const {
00438   if (IsUsable()){
00439     if (FmtStr==NULL){
00440       return TFlt::GetStr(GetByNm(MomNm));
00441     } else {
00442       return TFlt::GetStr(GetByNm(MomNm), FmtStr);
00443     }
00444   } else {
00445     return "X";
00446   }
00447 }
00448 
00449 TStr TMom::GetStr(
00450  const char& SepCh, const char& DelimCh,
00451  const bool& DecileP, const bool& PercentileP, const TStr& FmtStr) const {
00452   TChA ChA;
00453   if (IsUsable()){
00454     ChA+="["; ChA+=SepCh;
00455     ChA+="Vals"; ChA+=DelimCh; ChA+=TInt::GetStr(GetVals()); ChA+=SepCh;
00456     ChA+="Min"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMn(), FmtStr.CStr()); ChA+=SepCh;
00457     ChA+="Max"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMx(), FmtStr.CStr()); ChA+=SepCh;
00458     ChA+="Mean"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMean(), FmtStr.CStr()); ChA+=SepCh;
00459     //ChA+="Vari"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetVari(), FmtStr.CStr()); ChA+=SepCh;
00460     ChA+="SDev"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetSDev(), FmtStr.CStr()); ChA+=SepCh;
00461     //ChA+="SErr"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetSErr(), FmtStr.CStr()); ChA+=SepCh;
00462     ChA+="Quart1"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetQuart1(), FmtStr.CStr()); ChA+=SepCh;
00463     ChA+="Median"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetMedian(), FmtStr.CStr()); ChA+=SepCh;
00464     ChA+="Quart3"; ChA+=DelimCh; ChA+=TFlt::GetStr(GetQuart3(), FmtStr.CStr()); ChA+=SepCh;
00465     if (DecileP){
00466       for (int DecileN=0; DecileN<=10; DecileN++){
00467         ChA+="Dec"; ChA+=TInt::GetStr(DecileN);
00468         ChA+=DelimCh; ChA+=TFlt::GetStr(GetDecile(DecileN), FmtStr.CStr());
00469         ChA+=SepCh;
00470       }
00471     }
00472     if (PercentileP){
00473       for (int PercentileN=0; PercentileN<=100; PercentileN++){
00474         ChA+="Per"; ChA+=TInt::GetStr(PercentileN);
00475         ChA+=DelimCh; ChA+=TFlt::GetStr(GetPercentile(PercentileN), FmtStr.CStr());
00476         ChA+=SepCh;
00477       }
00478     }
00479     ChA+="]";
00480   } else {
00481     ChA="[Unusable]";
00482   }
00483   return ChA;
00484 }
00485 
00486 TStr TMom::GetNmVStr(const TStr& VarPfx,
00487  const char& SepCh, const bool& DecileP, const bool& PercentileP){
00488   TChA ChA;
00489   ChA+=VarPfx; ChA+="Vals"; ChA+=SepCh;
00490   ChA+=VarPfx; ChA+="Min"; ChA+=SepCh;
00491   ChA+=VarPfx; ChA+="Max"; ChA+=SepCh;
00492   ChA+=VarPfx; ChA+="Mean"; ChA+=SepCh;
00493   //ChA+=VarPfx; ChA+="Vari"; ChA+=SepCh;
00494   ChA+=VarPfx; ChA+="SDev"; ChA+=SepCh;
00495   //ChA+=VarPfx; ChA+="SErr"; ChA+=SepCh;
00496   ChA+=VarPfx; ChA+="Quart1"; ChA+=SepCh;
00497   ChA+=VarPfx; ChA+="Median"; ChA+=SepCh;
00498   ChA+=VarPfx; ChA+="Quart3";
00499   if (DecileP){
00500     ChA+=SepCh;
00501     for (int DecileN=0; DecileN<=10; DecileN++){
00502       ChA+=VarPfx; ChA+="Dec"; ChA+=TInt::GetStr(DecileN);
00503       if (DecileN<10){ChA+=SepCh;}
00504     }
00505   }
00506   if (PercentileP){
00507     ChA+=SepCh;
00508     for (int PercentileN=0; PercentileN<=100; PercentileN++){
00509       ChA+=VarPfx; ChA+="Per"; ChA+=TInt::GetStr(PercentileN);
00510       if (PercentileN<100){ChA+=SepCh;}
00511     }
00512   }
00513   return ChA;
00514 }
00515 
00516 TStr TMom::GetValVStr(
00517  const char& SepCh, const bool& DecileP, const bool& PercentileP) const {
00518   TChA ChA;
00519   if (IsUsable()){
00520     ChA+=TInt::GetStr(GetVals()); ChA+=SepCh;
00521     ChA+=TFlt::GetStr(GetMn()); ChA+=SepCh;
00522     ChA+=TFlt::GetStr(GetMx()); ChA+=SepCh;
00523     ChA+=TFlt::GetStr(GetMean()); ChA+=SepCh;
00524     //ChA+=TFlt::GetStr(GetVari()); ChA+=SepCh;
00525     ChA+=TFlt::GetStr(GetSDev()); ChA+=SepCh;
00526     //ChA+=TFlt::GetStr(GetSErr()); ChA+=SepCh;
00527     ChA+=TFlt::GetStr(GetQuart1()); ChA+=SepCh;
00528     ChA+=TFlt::GetStr(GetMedian()); ChA+=SepCh;
00529     ChA+=TFlt::GetStr(GetQuart3()); ChA+=SepCh;
00530     if (DecileP){
00531       for (int DecileN=0; DecileN<=10; DecileN++){
00532         ChA+=TFlt::GetStr(GetDecile(DecileN)); ChA+=SepCh;
00533       }
00534     }
00535     if (PercentileP){
00536       for (int PercentileN=0; PercentileN<=100; PercentileN++){
00537         ChA+=TFlt::GetStr(GetPercentile(PercentileN)); ChA+=SepCh;
00538       }
00539     }
00540   } else {
00541     int Vals=8;
00542     if (DecileP){Vals+=11;}
00543     if (PercentileP){Vals+=101;}
00544     for (int ValN=0; ValN<Vals; ValN++){
00545       ChA="[Unusable]";
00546       if (ValN<Vals-1){ChA+=SepCh;}
00547     }
00548   }
00549   return ChA;
00550 }
00551 
00553 // Correlation
00554 TCorr::TCorr(const TFltV& ValV1, const TFltV& ValV2):
00555   ValVLen(ValV1.Len()), CorrCf(), CorrCfPrb(), FisherZ(){
00556   static const double TINY=1.0e-20;
00557   IAssert(ValV1.Len()==ValV2.Len());
00558 
00559   // calculate the means
00560   double MeanVal1=0; double MeanVal2=0;
00561   {for (int ValN=0; ValN<ValVLen; ValN++){
00562     MeanVal1+=ValV1[ValN]; MeanVal2+=ValV2[ValN];}}
00563   MeanVal1/=ValVLen; MeanVal2/=ValVLen;
00564 
00565   // calculate correlation coefficient
00566   double yt, xt;
00567   double syy=0.0; double sxy=0.0; double sxx=0.0;
00568   {for (int ValN=0; ValN<ValVLen; ValN++){
00569     xt=ValV1[ValN]-MeanVal1;
00570     yt=ValV2[ValN]-MeanVal2;
00571     sxx+=xt*xt;
00572     syy+=yt*yt;
00573     sxy+=xt*yt;
00574   }}
00575   if (sxx*syy==0){
00576     CorrCf=0; //** not in numerical recipes - check why (pojavi se, ko so same nicle)
00577   } else {
00578     CorrCf=sxy/sqrt(sxx*syy);
00579   }
00580   // calculate correlation coefficient significance level
00581   double df=ValVLen-2;
00582   double t=CorrCf*sqrt(df/((1.0-CorrCf+TINY)*(1.0+CorrCf+TINY)));
00583   CorrCfPrb=TSpecFunc::BetaI(0.5*df,0.5,df/(df+t*t));
00584   // calculate Fisher's Z transformation
00585   FisherZ=0.5*log((1.0+(CorrCf)+TINY)/(1.0-(CorrCf)+TINY));
00586 }
00587 
00589 // Statistical Tests
00590 void TStatTest::AveVar(const TFltV& ValV, double& Ave, double& Var){
00591   Ave=0;
00592   for (int ValN=0; ValN<ValV.Len(); ValN++){
00593     Ave+=ValV[ValN];}
00594   Ave/=ValV.Len();
00595   Var=0;
00596   double ep=0;
00597   for (int ValN=0; ValN<ValV.Len(); ValN++){
00598     double s=ValV[ValN]-Ave;
00599     ep+=s;
00600     Var+=s*s;
00601   }
00602   Var=(Var-ep*ep/ValV.Len())/(ValV.Len()-1);
00603 }
00604 
00605 double TStatTest::KsProb(const double& Alam) {
00606   const double EPS1 = 0.001;
00607   const double EPS2 = 1.0e-8;
00608   double a2 = -2.0*Alam*Alam, fac = 2.0, sum = 0.0, term, termbf = 0.0;
00609   for (int j=1; j <= 100; j++) {
00610     term = fac*exp(a2*j*j);
00611     sum += term;
00612     if (fabs(term) <= EPS1*termbf || fabs(term) <= EPS2*sum)
00613       return sum;
00614     fac = -fac;
00615     termbf = fabs(term);
00616   }
00617   return 1.0;
00618 }
00619 
00620 void TStatTest::ChiSquareOne(
00621  const TFltV& ObservedBinV, const TFltV& ExpectedBinV,
00622  double& ChiSquareVal, double& SignificancePrb){
00623   IAssert(ObservedBinV.Len()==ExpectedBinV.Len());
00624   int Bins=ObservedBinV.Len();
00625   int Constraints=0;
00626   int DegreesOfFreedom=Bins-Constraints;
00627   ChiSquareVal=0.0;
00628   for (int BinN=0; BinN<Bins; BinN++){
00629     IAssert(ExpectedBinV[BinN]>0);
00630     double BinDiff=ObservedBinV[BinN]-ExpectedBinV[BinN];
00631     ChiSquareVal+=BinDiff*BinDiff/ExpectedBinV[BinN];
00632   }
00633   SignificancePrb=
00634    TSpecFunc::GammaQ(0.5*(DegreesOfFreedom), 0.5*(ChiSquareVal));
00635 }
00636 
00637 void TStatTest::ChiSquareTwo(
00638  const TFltV& ObservedBin1V, const TFltV& ObservedBin2V,
00639  double& ChiSquareVal, double& SignificancePrb){
00640   IAssert(ObservedBin1V.Len()==ObservedBin1V.Len());
00641   int Bins=ObservedBin1V.Len();
00642   int Constraints=0;
00643   int DegreesOfFreedom=Bins-Constraints;
00644   ChiSquareVal=0.0;
00645   for (int BinN=0; BinN<Bins; BinN++){
00646     if ((ObservedBin1V[BinN]==0.0) && (ObservedBin2V[BinN]==0.0)){
00647       DegreesOfFreedom--;
00648     } else {
00649       double BinDiff=ObservedBin1V[BinN]-ObservedBin2V[BinN];
00650       ChiSquareVal+=BinDiff*BinDiff/(ObservedBin1V[BinN]+ObservedBin2V[BinN]);
00651     }
00652   }
00653   SignificancePrb=
00654    TSpecFunc::GammaQ(0.5*(DegreesOfFreedom),0.5*(ChiSquareVal));
00655 }
00656 
00657 void TStatTest::TTest(
00658  const TFltV& ValV1, const TFltV& ValV2, double& TTestVal, double& TTestPrb){
00659   /*double Ave1; double Var1;
00660   AveVar(ValV1, Ave1, Var1);
00661   double Ave2; double Var2;
00662   AveVar(ValV2, Ave2, Var2);*/
00663 
00664   PMom Val1Mom=TMom::New(ValV1);
00665   PMom Val2Mom=TMom::New(ValV2);
00666   double ave1=Val1Mom->GetMean();
00667   double ave2=Val2Mom->GetMean();
00668   double var1=Val1Mom->GetVari();
00669   double var2=Val2Mom->GetVari();
00670   int n1=ValV1.Len();
00671   int n2=ValV2.Len();
00672 
00673   TTestVal=(ave1-ave2)/sqrt(var1/n1+var2/n2);
00674   double df=TMath::Sqr(var1/n1+var2/n2)/(TMath::Sqr(var1/n1)/(n1-1)+TMath::Sqr(var2/n2)/(n2-1));
00675   TTestPrb=TSpecFunc::BetaI(0.5*df, 0.5, df/(df+TMath::Sqr(TTestVal)));
00676 }
00677 
00678 void TStatTest::KsTest(const TFltV& ValV1, const TFltV& ValV2, double& DStat, double& PVal) {
00679   IAssert(ValV1.IsSorted() && ValV2.IsSorted());
00680   int i1=0, i2=0;
00681   double CumSum1=0.0, CumSum2=0.0, Cdf1=0.0, Cdf2=0.0;
00682   const double N1 = ValV1.Len();
00683   const double N2 = ValV2.Len();
00684   if (! (N1 > 0.0 && N2 > 0.0)) { DStat = 1.0;  PVal = 0.0;  return; }
00685   DStat=0.0; PVal=0.0;
00686   while (i1 < ValV1.Len() && i2 < ValV2.Len()) {
00687     const double X1 = ValV1[i1];
00688     const double X2 = ValV2[i2];
00689     if (X1 <= X2) {
00690       CumSum1 += 1;
00691       Cdf1 = (CumSum1 / N1);
00692       i1++;
00693     }
00694     if (X2 <= X1) {
00695       CumSum2 += 1;
00696       Cdf2 = (CumSum2 / N2);
00697       i2++;
00698     }
00699     DStat = TMath::Mx(DStat, fabs(Cdf1 - Cdf2));
00700   }
00701   const double En = sqrt( N1*N2 / (N1+N2));
00702   PVal = TStatTest::KsProb((En+0.12+0.11/En)*DStat);
00703 }
00704 
00705 void TStatTest::KsTest(const TFltPrV& ValCntV1, const TFltPrV& ValCntV2, double& DStat, double& PVal) {
00706   IAssert(ValCntV1.IsSorted() && ValCntV2.IsSorted());
00707   int i1=0, i2=0;
00708   double N1=0.0, N2=0.0, CumSum1=0.0, CumSum2=0.0, Cdf1=0.0, Cdf2=0.0;
00709   DStat=0.0;  PVal=0.0;
00710   for (int i = 0; i < ValCntV1.Len(); i++) N1 += ValCntV1[i].Val2;
00711   for (int i = 0; i < ValCntV2.Len(); i++) N2 += ValCntV2[i].Val2;
00712   if (! (N1 > 0.0 && N2 > 0.0)) { DStat = 1.0;  PVal = 0.0;  return; }
00713 
00714   while (i1 < ValCntV1.Len() && i2 < ValCntV2.Len()) {
00715     const double X1 = ValCntV1[i1].Val1;
00716     const double X2 = ValCntV2[i2].Val1;
00717     if (X1 <= X2) {
00718       CumSum1 += ValCntV1[i1].Val2;
00719       Cdf1 = (CumSum1 / N1);
00720       i1++;
00721     }
00722     if (X2 <= X1) {
00723       CumSum2 += ValCntV2[i2].Val2;
00724       Cdf2 = (CumSum2 / N2);
00725       i2++;
00726     }
00727     DStat = TMath::Mx(DStat, fabs(Cdf1 - Cdf2));
00728   }
00729   const double En = sqrt( N1*N2 / (N1+N2));
00730   PVal = TStatTest::KsProb((En+0.12+0.11/En)*DStat);
00731 }
00732 
00734 // Combinations
00735 bool TComb::GetNext(){
00736   if (ItemV.Len()==0){
00737     ItemV.Gen(Order, Order);
00738     for (int OrderN=0; OrderN<Order; OrderN++){
00739       ItemV[OrderN]=OrderN;}
00740     return true;
00741   } else {
00742     if (ItemV.Last()==Items-1){
00743       int OrderN=Order-1;
00744       while ((OrderN>=0)&&(ItemV[OrderN]==Items-(Order-OrderN-1)-1)){OrderN--;}
00745       if (OrderN<0){
00746         return false;
00747       } else {
00748         ItemV[OrderN]++;
00749         for (int SubOrderN=OrderN+1; SubOrderN<Order; SubOrderN++){
00750           ItemV[SubOrderN]=ItemV[SubOrderN-1]+1;}
00751         CombN++; return true;
00752       }
00753     } else {
00754       ItemV.Last()++; CombN++; return true;
00755     }
00756   }
00757 }
00758 
00759 int TComb::GetCombs() const {
00760   int LCombs=1; int HCombs=1;
00761   for (int OrderN=0; OrderN<Order; OrderN++){
00762     LCombs*=OrderN+1; HCombs*=Items-OrderN;}
00763   int Combs=HCombs/LCombs;
00764   return Combs;
00765 }
00766 
00767 void TComb::Wr(){
00768   printf("%d:[", GetCombN());
00769   for (int OrderN=0; OrderN<Order; OrderN++){
00770     if (OrderN>0){printf(" ");}
00771     printf("%d", ItemV[OrderN]());
00772   }
00773   printf("]\n");
00774 }
00775 
00777 // Linear-Regression
00778 PLinReg TLinReg::New(const TFltVV& _XVV, const TFltV& _YV, const TFltV& _SigV){
00779   PLinReg LinReg=PLinReg(new TLinReg());
00780   LinReg->XVV=_XVV;
00781   LinReg->YV=_YV;
00782   if (_SigV.Empty()){
00783     LinReg->SigV.Gen(LinReg->YV.Len());
00784     LinReg->SigV.PutAll(1);
00785   } else {
00786     LinReg->SigV=_SigV;
00787   }
00788   LinReg->Recs=LinReg->XVV.GetXDim();
00789   LinReg->Vars=LinReg->XVV.GetYDim();
00790   IAssert(LinReg->Recs>0);
00791   IAssert(LinReg->Vars>0);
00792   IAssert(LinReg->YV.Len()==LinReg->Recs);
00793   IAssert(LinReg->SigV.Len()==LinReg->Recs);
00794   LinReg->CovarVV.Gen(LinReg->Vars+1, LinReg->Vars+1);
00795   LinReg->CfV.Gen(LinReg->Vars+1);
00796   LinReg->NR_lfit();
00797   return LinReg;
00798 }
00799 
00800 void TLinReg::NR_covsrt(
00801  TFltVV& CovarVV, const int& Vars, const TIntV& ia, const int& mfit){
00802   for (int i=mfit+1; i<=Vars; i++){
00803     for (int j=1; j<=i; j++){
00804       CovarVV.At(i, j)=0; CovarVV.At(j, i)=0.0;}
00805   }
00806   int k=mfit;
00807   for (int j=Vars; j>=1; j--){
00808     if (ia[j]!=0){
00809       for (int i=1; i<=Vars; i++){Swap(CovarVV.At(i, k), CovarVV.At(i, j));}
00810       {for (int i=1; i<=Vars; i++){Swap(CovarVV.At(k, i), CovarVV.At(j, i));}}
00811       k--;
00812     }
00813   }
00814 }
00815 
00816 void TLinReg::NR_gaussj(TFltVV& a, const int& n, TFltVV& b, const int& m){
00817   int i, icol=0, irow=0, j, k, l, ll;
00818   double big, dum, pivinv;
00819 
00820   TIntV indxc(n+1);
00821   TIntV indxr(n+1);
00822   TIntV ipiv(n+1);
00823   for (j=1; j<=n; j++){ipiv[j]=0;}
00824   for (i=1; i<=n; i++){
00825     big=0.0;
00826     for (j=1; j<=n; j++){
00827       if (ipiv[j]!=1){
00828         for (k=1; k<=n; k++){
00829           if (ipiv[k]==0){
00830             if (fabs(double(a.At(j, k))) >= big){
00831               big=fabs(double(a.At(j, k)));
00832               irow=j;
00833               icol=k;
00834             }
00835           } else
00836           if (ipiv[k]>1){
00837             TExcept::Throw("Singular Matrix(1) in Gauss");}
00838         }
00839       }
00840     }
00841     ipiv[icol]++;
00842     if (irow != icol){
00843       for (l=1; l<=n; l++){Swap(a.At(irow, l), a.At(icol, l));}
00844       for (l=1; l<=m; l++){Swap(b.At(irow, l), b.At(icol, l));}
00845     }
00846     indxr[i]=irow;
00847     indxc[i]=icol;
00848     if (a.At(icol, icol)==0.0){
00849       TExcept::Throw("Singular Matrix(1) in Gauss");}
00850     pivinv=1.0/a.At(icol, icol);
00851     a.At(icol, icol)=1.0;
00852     for (l=1; l<=n; l++){a.At(icol, l)=a.At(icol, l)*pivinv;}
00853     for (l=1; l<=m; l++){b.At(icol, l)=b.At(icol, l)*pivinv;}
00854     for (ll=1; ll<=n; ll++){
00855       if (ll != icol){
00856         dum=a.At(ll, icol);
00857         a.At(ll, icol)=0.0;
00858         for (l=1;l<=n;l++){a.At(ll, l)-=a.At(icol, l)*dum;}
00859         for (l=1;l<=m;l++){b.At(ll, l)-=b.At(icol, l)*dum;}
00860       }
00861     }
00862   }
00863   for (l=n; l>=1; l--){
00864     if (indxr[l]!=indxc[l]){
00865       for (k=1; k<=n; k++){
00866         Swap(a.At(k, indxr[l]), a.At(k, indxc[l]));}
00867     }
00868   }
00869 }
00870 
00871 void TLinReg::NR_lfit(){
00872   int i,j,k,l,m,mfit=0;
00873   double ym,wt,sum,sig2i;
00874 
00875   TIntV ia(Vars+1); for (i=1; i<=Vars; i++){ia[i]=1;}
00876   TFltVV beta(Vars+1, 1+1);
00877   TFltV afunc(Vars+1);
00878   for (j=1;j<=Vars;j++){
00879     if (ia[j]!=0){mfit++;}}
00880   if (mfit==0){TExcept::Throw("No parameters to be fitted in LFit");}
00881   for (j=1; j<=mfit; j++){
00882     for (k=1; k<=mfit; k++){CovarVV.At(j, k)=0.0;}
00883     beta.At(j, 1)=0.0;
00884   }
00885   for (i=1; i<=Recs; i++){
00886     GetXV(i, afunc); // funcs(XVV[i],afunc,Vars);
00887     ym=GetY(i);
00888     if (mfit<Vars){
00889       for (j=1;j<=Vars;j++){
00890         if (ia[j]==0){ym-=CfV[j]*afunc[j];}}
00891     }
00892     sig2i=1.0/TMath::Sqr(GetSig(i));
00893     for (j=0, l=1; l<=Vars; l++){
00894       if (ia[l]!=0){
00895         wt=afunc[l]*sig2i;
00896         for (j++, k=0, m=1; m<=l; m++){
00897           if (ia[m]!=0){CovarVV.At(j, ++k)+=wt*afunc[m];}
00898         }
00899         beta.At(j, 1)+=ym*wt;
00900       }
00901     }
00902   }
00903   for (j=2; j<=mfit; j++){
00904     for (k=1; k<j; k++){CovarVV.At(k, j)=CovarVV.At(j, k);}
00905   }
00906   NR_gaussj(CovarVV, mfit, beta, 1);
00907   for (j=0, l=1; l<=Vars; l++){
00908     if (ia[l]!=0){CfV[l]=beta.At(++j, 1);}
00909   }
00910   ChiSq=0.0;
00911   for (i=1; i<=Recs; i++){
00912     GetXV(i, afunc); // funcs(XVV[i],afunc,Vars);
00913     for (sum=0.0, j=1; j<=Vars; j++){sum+=CfV[j]*afunc[j];}
00914     ChiSq+=TMath::Sqr((GetY(i)-sum)/GetSig(i));
00915   }
00916   NR_covsrt(CovarVV, Vars, ia, mfit);
00917 }
00918 
00919 void TLinReg::Wr() const {
00920   printf("\n%11s %21s\n","parameter","uncertainty");
00921   for (int i=0; i<Vars;i++){
00922     printf("  a[%1d] = %8.6f %12.6f\n",
00923      i+1, GetCf(i), GetCfUncer(i));
00924   }
00925   printf("chi-squared = %12f\n", GetChiSq());
00926 
00927   printf("full covariance matrix\n");
00928   {for (int i=0;i<Vars;i++){
00929     for (int j=0;j<Vars;j++){
00930       printf("%12f", GetCovar(i, j));}
00931     printf("\n");
00932   }}
00933 }
00934 
00936 // Singular-Value-Decomposition
00937 PSvd TSvd::New(const TFltVV& _XVV, const TFltV& _YV, const TFltV& _SigV){
00938   PSvd Svd=PSvd(new TSvd());
00939   Svd->XVV=_XVV;
00940   Svd->YV=_YV;
00941   if (_SigV.Empty()){
00942     Svd->SigV.Gen(Svd->YV.Len());
00943     Svd->SigV.PutAll(1);
00944   } else {
00945     Svd->SigV=_SigV;
00946   }
00947   Svd->Recs=Svd->XVV.GetXDim();
00948   Svd->Vars=Svd->XVV.GetYDim();
00949   IAssert(Svd->Recs>0);
00950   IAssert(Svd->Vars>0);
00951   IAssert(Svd->YV.Len()==Svd->Recs);
00952   IAssert(Svd->SigV.Len()==Svd->Recs);
00953   Svd->CovarVV.Gen(Svd->Vars+1, Svd->Vars+1);
00954   Svd->CfV.Gen(Svd->Vars+1);
00955   Svd->NR_svdfit();
00956   return Svd;
00957 }
00958 
00959 double TSvd::NR_pythag(double a, double b){
00960   double absa,absb;
00961   absa=fabs(a);
00962   absb=fabs(b);
00963   if (absa > absb){
00964     return absa*sqrt(1.0+TMath::Sqr(absb/absa));
00965   } else {
00966     return (absb == 0.0 ? 0.0 : absb*sqrt(1.0+TMath::Sqr(absa/absb)));
00967   }
00968 }
00969 
00970 void TSvd::NR_svdcmp(TFltVV& a, int m, int n, TFltV& w, TFltVV& v){
00971   int flag,i,its,j,jj,k,l=0,nm;
00972   double anorm,c,f,g,h,s,scale,x,y,z;
00973 
00974   TFltV rv1(n+1);
00975   g=scale=anorm=0.0;
00976   for (i=1;i<=n;i++) {
00977     l=i+1;
00978     rv1[i]=scale*g;
00979     g=s=scale=0.0;
00980     if (i <= m) {
00981       for (k=i;k<=m;k++) scale += fabs(double(a.At(k,i)));
00982       if (scale) {
00983         for (k=i;k<=m;k++) {
00984           a.At(k,i) /= scale;
00985           s += a.At(k,i)*a.At(k,i);
00986         }
00987         f=a.At(i,i);
00988         g = -NR_SIGN(sqrt(s),f);
00989         h=f*g-s;
00990         a.At(i,i)=f-g;
00991         for (j=l;j<=n;j++) {
00992           for (s=0.0,k=i;k<=m;k++) s += a.At(k,i)*a(k,j);
00993           f=s/h;
00994           for (k=i;k<=m;k++) a.At(k,j) += f*a.At(k,i);
00995         }
00996         for (k=i;k<=m;k++) a.At(k,i) *= scale;
00997       }
00998     }
00999     w[i]=scale *g;
01000     g=s=scale=0.0;
01001     if (i <= m && i != n) {
01002       for (k=l;k<=n;k++) scale += fabs(double(a.At(i,k)));
01003       if (scale) {
01004         for (k=l;k<=n;k++) {
01005           a.At(i,k) /= scale;
01006           s += a.At(i,k)*a.At(i,k);
01007         }
01008         f=a.At(i,l);
01009         g = -NR_SIGN(sqrt(s),f);
01010         h=f*g-s;
01011         a.At(i,l)=f-g;
01012         for (k=l;k<=n;k++) rv1[k]=a.At(i,k)/h;
01013         for (j=l;j<=m;j++) {
01014           for (s=0.0,k=l;k<=n;k++) s += a.At(j,k)*a.At(i,k);
01015           for (k=l;k<=n;k++) a.At(j,k) += s*rv1[k];
01016         }
01017         for (k=l;k<=n;k++) a.At(i,k) *= scale;
01018       }
01019     }
01020     anorm=NR_FMAX(anorm,(fabs(double(w[i]))+fabs(double(rv1[i]))));
01021   }
01022   for (i=n;i>=1;i--) {
01023     if (i < n) {
01024       if (g) {
01025         for (j=l;j<=n;j++)
01026           v.At(j,i)=(a.At(i,j)/a.At(i,l))/g;
01027         for (j=l;j<=n;j++) {
01028           for (s=0.0,k=l;k<=n;k++) s += a.At(i,k)*v.At(k,j);
01029           for (k=l;k<=n;k++) v.At(k,j) += s*v.At(k,i);
01030         }
01031       }
01032       for (j=l;j<=n;j++) v.At(i,j)=v.At(j,i)=0.0;
01033     }
01034     v.At(i,i)=1.0;
01035     g=rv1[i];
01036     l=i;
01037   }
01038   for (i=NR_IMIN(m,n);i>=1;i--) {
01039     l=i+1;
01040     g=w[i];
01041     for (j=l;j<=n;j++) a.At(i,j)=0.0;
01042     if (g) {
01043       g=1.0/g;
01044       for (j=l;j<=n;j++) {
01045         for (s=0.0,k=l;k<=m;k++) s += a.At(k,i)*a.At(k,j);
01046         f=(s/a.At(i,i))*g;
01047         for (k=i;k<=m;k++) a.At(k,j) += f*a.At(k,i);
01048       }
01049       for (j=i;j<=m;j++) a.At(j,i) *= g;
01050     } else for (j=i;j<=m;j++) a.At(j,i)=0.0;
01051     a.At(i,i)++;
01052   }
01053   for (k=n;k>=1;k--) {
01054     for (its=1;its<=30;its++) {
01055       flag=1;
01056       for (l=k;l>=1;l--) {
01057         nm=l-1;
01058         if ((double)(fabs(double(rv1[l])+anorm)) == anorm) {
01059           flag=0;
01060           break;
01061         }
01062         if ((double)(fabs(double(w[nm]))+anorm) == anorm) break;
01063       }
01064       if (flag) {
01065         c=0.0;
01066         s=1.0;
01067         for (i=l;i<=k;i++) {
01068           f=s*rv1[i];
01069           rv1[i]=c*rv1[i];
01070           if ((double)(fabs(f)+anorm) == anorm) break;
01071           g=w[i];
01072           h=NR_pythag(f,g);
01073           w[i]=h;
01074           h=1.0/h;
01075           c=g*h;
01076           s = -f*h;
01077           for (j=1;j<=m;j++) {
01078             y=a.At(j,nm);
01079             z=a.At(j,i);
01080             a.At(j,nm)=y*c+z*s;
01081             a.At(j,i)=z*c-y*s;
01082           }
01083         }
01084       }
01085       z=w[k];
01086       if (l == k) {
01087         if (z < 0.0) {
01088           w[k] = -z;
01089           for (j=1;j<=n;j++) v.At(j,k) = -v.At(j,k);
01090         }
01091         break;
01092       }
01093       if (its==30){
01094         TExcept::Throw("no convergence in 30 svdcmp iterations");}
01095       x=w[l];
01096       nm=k-1;
01097       y=w[nm];
01098       g=rv1[nm];
01099       h=rv1[k];
01100       f=((y-z)*(y+z)+(g-h)*(g+h))/(2.0*h*y);
01101       g=NR_pythag(f,1.0);
01102       f=((x-z)*(x+z)+h*((y/(f+NR_SIGN(g,f)))-h))/x;
01103       c=s=1.0;
01104       for (j=l;j<=nm;j++) {
01105         i=j+1;
01106         g=rv1[i];
01107         y=w[i];
01108         h=s*g;
01109         g=c*g;
01110         z=NR_pythag(f,h);
01111         rv1[j]=z;
01112         c=f/z;
01113         s=h/z;
01114         f=x*c+g*s;
01115         g = g*c-x*s;
01116         h=y*s;
01117         y *= c;
01118         for (jj=1;jj<=n;jj++) {
01119           x=v.At(jj,j);
01120           z=v.At(jj,i);
01121           v.At(jj,j)=x*c+z*s;
01122           v.At(jj,i)=z*c-x*s;
01123         }
01124         z=NR_pythag(f,h);
01125         w[j]=z;
01126         if (z) {
01127           z=1.0/z;
01128           c=f*z;
01129           s=h*z;
01130         }
01131         f=c*g+s*y;
01132         x=c*y-s*g;
01133         for (jj=1;jj<=m;jj++) {
01134           y=a.At(jj,j);
01135           z=a.At(jj,i);
01136           a.At(jj,j)=y*c+z*s;
01137           a.At(jj,i)=z*c-y*s;
01138         }
01139       }
01140       rv1[l]=0.0;
01141       rv1[k]=f;
01142       w[k]=x;
01143     }
01144   }
01145 }
01146 
01147 void TSvd::NR_svbksb(
01148  TFltVV& u, TFltV& w, TFltVV& v, int m, int n, TFltV& b, TFltV& x){
01149   int jj,j,i;
01150   double s;
01151 
01152   TFltV tmp(n+1);
01153   for (j=1;j<=n;j++) {
01154     s=0.0;
01155     if (w[j]) {
01156       for (i=1;i<=m;i++) s += u.At(i,j)*b[i];
01157       s /= w[j];
01158     }
01159     tmp[j]=s;
01160   }
01161   for (j=1;j<=n;j++) {
01162     s=0.0;
01163     for (jj=1;jj<=n;jj++) s += v.At(j,jj)*tmp[jj];
01164     x[j]=s;
01165   }
01166 }
01167 
01168 void TSvd::NR_svdvar(TFltVV& v, int ma, TFltV& w, TFltVV& cvm){
01169   int k,j,i;
01170   double sum;
01171 
01172   TFltV wti(ma+1);
01173   for (i=1;i<=ma;i++) {
01174     wti[i]=0.0;
01175     if (w[i]) wti[i]=1.0/(w[i]*w[i]);
01176   }
01177   for (i=1;i<=ma;i++) {
01178     for (j=1;j<=i;j++) {
01179       for (sum=0.0,k=1;k<=ma;k++) sum += v.At(i,k)*v.At(j,k)*wti[k];
01180       cvm.At(j,i)=cvm.At(i,j)=sum;
01181     }
01182   }
01183 }
01184 
01185 void TSvd::NR_svdfit(){
01186   int j,i;
01187   double wmax,tmp,thresh,sum;
01188   double TOL=1.0e-5;
01189 
01190   TFltVV u(Recs+1, Vars+1);
01191   TFltVV v(Vars+1, Vars+1);
01192   TFltV w(Vars+1);
01193   TFltV b(Recs+1);
01194   TFltV afunc(Vars+1);
01195   for (i=1;i<=Recs;i++) {
01196     GetXV(i, afunc); // (*funcs)(x[i],afunc,Vars);
01197     tmp=1.0/GetSig(i);
01198     for (j=1;j<=Vars;j++){u.At(i,j)=afunc[j]*tmp;}
01199     b[i]=GetY(i)*tmp;
01200   }
01201   NR_svdcmp(u,Recs,Vars,w,v);
01202   wmax=0.0;
01203   for (j=1;j<=Vars;j++){
01204     if (w[j] > wmax){wmax=w[j];}}
01205   thresh=TOL*wmax;
01206   for (j=1;j<=Vars;j++){
01207     if (double(w[j])<thresh){w[j]=0.0;}}
01208   NR_svbksb(u,w,v,Recs,Vars,b,CfV);
01209   ChiSq=0.0;
01210   for (i=1;i<=Recs;i++) {
01211     GetXV(i, afunc); // (*funcs)(x[i],afunc,Vars);
01212     for (sum=0.0,j=1;j<=Vars;j++){sum += CfV[j]*afunc[j];}
01213     ChiSq += (tmp=(GetY(i)-sum)/GetSig(i),tmp*tmp);
01214   }
01215   // covariance matrix calculation
01216   CovarVV.Gen(Vars+1, Vars+1);
01217   NR_svdvar(v, Vars, w, CovarVV);
01218 }
01219 
01220 void TSvd::GetCfV(TFltV& _CfV){
01221   _CfV=CfV; _CfV.Del(0);
01222 }
01223 
01224 void TSvd::GetCfUncerV(TFltV& CfUncerV){
01225   CfUncerV.Gen(Vars);
01226   for (int VarN=0; VarN<Vars; VarN++){
01227     CfUncerV[VarN]=GetCfUncer(VarN);
01228   }
01229 }
01230 
01231 // all vectors (matrices) start with 0 
01232 void TSvd::Svd(const TFltVV& InMtx, TFltVV& LSingV, TFltV& SingValV, TFltVV& RSingV) {
01233   //LSingV = InMtx;
01234   LSingV.Gen(InMtx.GetYDim()+1, InMtx.GetYDim()+1);
01235   // create 1 based adjacency matrix
01236   for (int x = 0; x < InMtx.GetXDim(); x++) {
01237     for (int y = 0; y < InMtx.GetYDim(); y++) {
01238       LSingV.At(x+1, y+1) = InMtx.At(x, y);
01239     }
01240   }
01241   RSingV.Gen(InMtx.GetYDim()+1, InMtx.GetYDim()+1);
01242   SingValV.Gen(InMtx.GetYDim()+1);
01243   TSvd::NR_svdcmp(LSingV, InMtx.GetXDim(), InMtx.GetYDim(), SingValV, RSingV);
01244   // 0-th singular value/vector is full of zeros, delete it
01245   SingValV.Del(0);
01246   LSingV.DelX(0); LSingV.DelY(0);
01247   RSingV.DelX(0); RSingV.DelY(0);
01248 }
01249 
01250 // InMtx1 is 1-based (0-th row/column are full of zeros)
01251 // returned singular vectors/values are 0 based
01252 void TSvd::Svd1Based(const TFltVV& InMtx1, TFltVV& LSingV, TFltV& SingValV, TFltVV& RSingV) {
01253   LSingV = InMtx1;
01254   SingValV.Gen(InMtx1.GetYDim());
01255   RSingV.Gen(InMtx1.GetYDim(), InMtx1.GetYDim());
01256   TSvd::NR_svdcmp(LSingV, InMtx1.GetXDim()-1, InMtx1.GetYDim()-1, SingValV, RSingV);
01257   // 0-th singular value/vector is full of zeros, delete it
01258   SingValV.Del(0);
01259   LSingV.DelX(0); LSingV.DelY(0);
01260   RSingV.DelX(0); RSingV.DelY(0);
01261 }
01262 
01263 void TSvd::Wr() const {
01264   printf("\n%11s %21s\n","parameter","uncertainty");
01265   for (int i=0; i<Vars;i++){
01266     printf("  a[%1d] = %8.6f %12.6f\n",
01267      i+1, GetCf(i), GetCfUncer(i));
01268   }
01269   printf("chi-squared = %12f\n", GetChiSq());
01270 
01271   printf("full covariance matrix\n");
01272   {for (int i=0;i<Vars;i++){
01273     for (int j=0;j<Vars;j++){
01274       printf("%12f", GetCovar(i, j));}
01275     printf("\n");
01276   }}
01277 }
01278 
01280 // Histogram
01281 void THist::Add(const double& Val, const bool& OnlyInP) {
01282         // get bucket number
01283     const int BucketN = int(floor((Val - MnVal) / BucketSize));
01284         if (OnlyInP) { 
01285                 // value should be inside
01286                 EAssert(MnVal <= Val && Val <= MxVal);
01287                 BucketV[BucketN]++;
01288         } else {
01289                 // value does not need to be inside
01290                 if (BucketN < 0) {
01291                         BucketV[0]++;
01292                 } else if (BucketN < BucketV.Len()) {
01293                         BucketV[BucketN]++;
01294                 } else {
01295                         BucketV.Last()++;
01296                 }
01297         }
01298         // for computing percentage
01299         Vals++;
01300 }
01301 
01302 void THist::SaveStat(const TStr& ValNm, TSOut& FOut) const {
01303     FOut.PutStrLn("#" + ValNm + ": " + Vals.GetStr());
01304     const int Buckets = BucketV.Len() - 1;
01305     for (int BucketN = 0; BucketN < Buckets; BucketN++) {
01306         FOut.PutStrLn(TStr::Fmt("%d-%d\t%d", BucketSize*BucketN,
01307             BucketSize*(BucketN+1), BucketV[BucketN]()));
01308     }
01309     if (BucketV.Last() > 0) {
01310         FOut.PutStrLn(TStr::Fmt("%d-\t%d", BucketSize*Buckets, BucketV.Last()()));
01311     }
01312 
01313 }
01314