TH1K.cxx

Go to the documentation of this file.
00001 // @(#)root/hist:$Id: TH1K.cxx 35245 2010-09-13 14:03:43Z brun $
00002 // Author: Victor Perevoztchikov <perev@bnl.gov>  21/02/2001
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2000, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 #include <stdlib.h>
00013 
00014 #include "Riostream.h"
00015 #include "TH1K.h"
00016 #include "TMath.h"
00017 
00018 //______________________________________________________________________________
00019 // The TH1K Class
00020 //
00021 //  TH1K class supports the nearest K Neighbours method,
00022 //       widely used in cluster analysis.
00023 //       This method is especially useful for small statistics.
00024 //
00025 //       In this method :
00026 //         DensityOfProbability ~ 1/DistanceToNearestKthNeighbour
00027 //
00028 //      Ctr TH1K::TH1K(name,title,nbins,xlow,xup,K=0)
00029 //      differs from TH1F only by "K"
00030 //      K - is the order of K Neighbours method, usually >=3
00031 //      K = 0, means default, where K is selected by TH1K in such a way
00032 //             that DistanceToNearestKthNeighbour > BinWidth and K >=3
00033 //
00034 //  This class has been implemented by Victor Perevoztchikov <perev@bnl.gov>
00035 
00036 ClassImp(TH1K)
00037 
00038 
00039 //______________________________________________________________________________
00040 TH1K::TH1K(): TH1(), TArrayF()
00041 {
00042    // Constructor.
00043 
00044    fDimension = 1;
00045    fNIn   = 0;
00046    fReady = 0;
00047    fKOrd  = 3;
00048    fKCur  = 0;
00049 }
00050 
00051 
00052 //______________________________________________________________________________
00053 TH1K::TH1K(const char *name,const char *title,Int_t nbins,Double_t xlow,Double_t xup,Int_t k)
00054      : TH1(name,title,nbins,xlow,xup), TArrayF(100)
00055 {
00056    // Create a 1-Dim histogram with fix bins of type float
00057    // (see TH1K::TH1 for explanation of parameters)
00058 
00059    fDimension = 1;
00060    fNIn   = 0;
00061    fReady = 0;
00062    fKOrd  = k;
00063    fKCur  = 0;
00064 }
00065 
00066 
00067 //______________________________________________________________________________
00068 TH1K::~TH1K()
00069 {
00070    // Destructor.
00071 }
00072 
00073 
00074 //______________________________________________________________________________
00075 Int_t TH1K::Fill(Double_t x)
00076 {
00077    // Increment bin with abscissa X by 1.
00078    //
00079    // if x is less than the low-edge of the first bin, the Underflow bin is incremented
00080    // if x is greater than the upper edge of last bin, the Overflow bin is incremented
00081    //
00082    // If the storage of the sum of squares of weights has been triggered,
00083    // via the function Sumw2, then the sum of the squares of weights is incremented
00084    // by 1 in the bin corresponding to x.
00085 
00086    fReady = 0;
00087    Int_t bin;
00088    fEntries++;
00089    bin =fXaxis.FindBin(x);
00090    if (bin == 0 || bin > fXaxis.GetNbins()) {
00091       if (!fgStatOverflows) return -1;
00092    }
00093    ++fTsumw;
00094    ++fTsumw2;
00095    fTsumwx  += x;
00096    fTsumwx2 += x*x;
00097    fReady = 0;
00098    if (fNIn == fN) Set(fN*2);
00099    AddAt(x,fNIn++);
00100    return bin;
00101 }
00102 
00103 
00104 //______________________________________________________________________________
00105 Double_t TH1K::GetBinContent(Int_t bin) const
00106 {
00107    // Return content of global bin number bin.
00108 
00109    if (!fReady) {
00110       ((TH1K*)this)->Sort();
00111       ((TH1K*)this)->fReady=1;
00112    }
00113    if (!fNIn)   return 0.;
00114    float x = GetBinCenter(bin);
00115    int left = TMath::BinarySearch(fNIn,fArray,x);
00116    int jl=left,jr=left+1,nk,nkmax =fKOrd;
00117    float fl,fr,ff=0.,ffmin=1.e-6;
00118    if (!nkmax) {nkmax = 3; ffmin = GetBinWidth(bin);}
00119    if (nkmax >= fNIn) nkmax = fNIn-1;
00120    for (nk = 1; nk <= nkmax || ff <= ffmin; nk++) {
00121       fl = (jl>=0  ) ? TMath::Abs(fArray[jl]-x) : 1.e+20;
00122       fr = (jr<fNIn) ? TMath::Abs(fArray[jr]-x) : 1.e+20;
00123       if (jl<0 && jr>=fNIn) break;
00124       if (fl < fr) { ff = fl; jl--;}
00125       else         { ff = fr; jr++;}
00126    }
00127    ((TH1K*)this)->fKCur = nk - 1;
00128    return fNIn * 0.5*fKCur/((float)(fNIn+1))*GetBinWidth(bin)/ff;
00129 }
00130 
00131 
00132 //______________________________________________________________________________
00133 Double_t TH1K::GetBinError(Int_t bin) const
00134 {
00135    // Return content of global bin error.
00136 
00137    return TMath::Sqrt(((double)(fNIn-fKCur+1))/((fNIn+1)*(fKCur-1)))*GetBinContent(bin);
00138 }
00139 
00140 
00141 //______________________________________________________________________________
00142 void   TH1K::Reset(Option_t *option)
00143 {
00144    // Reset.
00145 
00146    fNIn   =0;
00147    fReady = 0;
00148    TH1::Reset(option);
00149 }
00150 
00151 
00152 //______________________________________________________________________________
00153 void TH1K::SavePrimitive(ostream &out, Option_t *option /*= ""*/)
00154 {
00155    // Save primitive as a C++ statement(s) on output stream out
00156    // Note the following restrictions in the code generated:
00157    //  - variable bin size not implemented
00158    //  - Objects in list of functions not saved (fits)
00159    //  - Contours not saved
00160 
00161    char quote = '"';
00162    out<<"   "<<endl;
00163    out<<"   "<<"TH1 *";
00164 
00165    out<<GetName()<<" = new "<<ClassName()<<"("<<quote<<GetName()<<quote<<","<<quote<<GetTitle()<<quote
00166                  <<","<<GetXaxis()->GetNbins()
00167                  <<","<<GetXaxis()->GetXmin()
00168                  <<","<<GetXaxis()->GetXmax()
00169                  <<","<<fKOrd;
00170    out <<");"<<endl;
00171 
00172    if (fDirectory == 0) {
00173       out<<"   "<<GetName()<<"->SetDirectory(0);"<<endl;
00174    }
00175    if (TestBit(kNoStats)) {
00176       out<<"   "<<GetName()<<"->SetStats(0);"<<endl;
00177    }
00178    if (fOption.Length() != 0) {
00179       out<<"   "<<GetName()<<"->SetOption("<<quote<<fOption.Data()<<quote<<");"<<endl;
00180    }
00181 
00182    if (fNIn) {
00183       out<<"   Float_t Arr[]={"<<endl;
00184       for (int i=0; i<fNIn; i++) {
00185          out<< fArray[i];
00186          if (i != fNIn-1) {out<< ",";} else { out<< "};";}
00187          if (i%10 == 9)   {out<< endl;}
00188       }
00189       out<< endl;
00190       out<<"   for(int i=0;i<" << fNIn << ";i++)"<<GetName()<<"->Fill(Arr[i]);";
00191       out<< endl;
00192    }
00193    SaveFillAttributes(out,GetName(),0,1001);
00194    SaveLineAttributes(out,GetName(),1,1,1);
00195    SaveMarkerAttributes(out,GetName(),1,1,1);
00196    fXaxis.SaveAttributes(out,GetName(),"->GetXaxis()");
00197    fYaxis.SaveAttributes(out,GetName(),"->GetYaxis()");
00198    fZaxis.SaveAttributes(out,GetName(),"->GetZaxis()");
00199    TString opt = option;
00200    opt.ToLower();
00201    if (!opt.Contains("nodraw")) {
00202       out<<"   "<<GetName()<<"->Draw("
00203       <<quote<<option<<quote<<");"<<endl;
00204    }
00205 }
00206 
00207 
00208 //______________________________________________________________________________
00209 static int TH1K_fcompare(const void *f1,const void *f2)
00210 {
00211    // Compare.
00212 
00213    if (*((float*)f1) < *((float*)f2)) return -1;
00214    if (*((float*)f1) > *((float*)f2)) return  1;
00215    return 0;
00216 }
00217 
00218 
00219 //______________________________________________________________________________
00220 void TH1K::Sort()
00221 {
00222    // Sort.
00223 
00224    if (fNIn<2) return;
00225    qsort(GetArray(),fNIn,sizeof(Float_t),&TH1K_fcompare);
00226 }

Generated on Tue Jul 5 14:24:19 2011 for ROOT_528-00b_version by  doxygen 1.5.1