TFoamMaxwt.cxx

Go to the documentation of this file.
00001 // @(#)root/foam:$Id: TFoamMaxwt.cxx 22726 2008-03-19 09:53:41Z pcanal $
00002 // Author: S. Jadach <mailto:Stanislaw.jadach@ifj.edu.pl>, P.Sawicki <mailto:Pawel.Sawicki@ifj.edu.pl>
00003 
00004 //____________________________________________________________________________
00005 //
00006 // Class  TFoamMaxwt
00007 // =================
00008 // Small auxiliary class for controlling MC weight.
00009 // It provides certain measure of the "maximum weight"
00010 // depending on small user-parameter "epsilon".
00011 // It creates and uses 2 histograms of the TH1D class.
00012 // User defines no. of bins nBin,  nBin=1000 is  recommended
00013 // wmax defines weight range (1,wmax), it is adjusted "manually"
00014 //
00015 //____________________________________________________________________________
00016 
00017 
00018 #include "Riostream.h"
00019 #include "TH1.h"
00020 #include "TFoamMaxwt.h"
00021 
00022 ClassImp(TFoamMaxwt);
00023 
00024 //____________________________________________________________________________
00025 TFoamMaxwt::TFoamMaxwt()
00026 {
00027 // Constructor for streamer
00028    fNent = 0;
00029    fnBin = 0;
00030    fWtHst1 = 0;
00031    fWtHst2 = 0;
00032 }
00033 
00034 //____________________________________________________________________________
00035 TFoamMaxwt::TFoamMaxwt(Double_t wmax, Int_t nBin)
00036 {
00037 // Principal user constructor
00038    fNent = 0;
00039    fnBin = nBin;
00040    fwmax = wmax;
00041    fWtHst1 = new TH1D("TFoamMaxwt_hst_Wt1","Histo of weight   ",nBin,0.0,wmax);
00042    fWtHst2 = new TH1D("TFoamMaxwt_hst_Wt2","Histo of weight**2",nBin,0.0,wmax);
00043    fWtHst1->SetDirectory(0);// exclude from diskfile
00044    fWtHst2->SetDirectory(0);// and enable deleting
00045 }
00046 
00047 //______________________________________________________________________________
00048 TFoamMaxwt::TFoamMaxwt(TFoamMaxwt &From): TObject(From)
00049 {
00050 // Explicit COPY CONSTRUCTOR (unused, so far)
00051    fnBin   = From.fnBin;
00052    fwmax   = From.fwmax;
00053    fWtHst1 = From.fWtHst1;
00054    fWtHst2 = From.fWtHst2;
00055    Error("TFoamMaxwt","COPY CONSTRUCTOR NOT TESTED!");
00056 }
00057 
00058 //_______________________________________________________________________________
00059 TFoamMaxwt::~TFoamMaxwt()
00060 {
00061 // Destructor
00062    delete fWtHst1; // For this SetDirectory(0) is needed!
00063    delete fWtHst2; //
00064    fWtHst1=0;
00065    fWtHst2=0;
00066 }
00067 //_______________________________________________________________________________
00068 void TFoamMaxwt::Reset()
00069 {
00070 // Reseting weight analysis
00071    fNent = 0;
00072    fWtHst1->Reset();
00073    fWtHst2->Reset();
00074 }
00075 
00076 //_______________________________________________________________________________
00077 TFoamMaxwt& TFoamMaxwt::operator=(const TFoamMaxwt &From)
00078 {
00079 // substitution =
00080    if (&From == this) return *this;
00081    fnBin = From.fnBin;
00082    fwmax = From.fwmax;
00083    fWtHst1 = From.fWtHst1;
00084    fWtHst2 = From.fWtHst2;
00085    return *this;
00086 }
00087 
00088 //________________________________________________________________________________
00089 void TFoamMaxwt::Fill(Double_t wt)
00090 {
00091 // Filling analyzed weight
00092    fNent =  fNent+1.0;
00093    fWtHst1->Fill(wt,1.0);
00094    fWtHst2->Fill(wt,wt);
00095 }
00096 
00097 //________________________________________________________________________________
00098 void TFoamMaxwt::Make(Double_t eps, Double_t &MCeff)
00099 {
00100 // Calculates Efficiency= aveWt/wtLim for a given tolerance level epsilon<<1
00101 // To be called at the end of the MC run.
00102 
00103    Double_t wtLim,aveWt;
00104    GetMCeff(eps, MCeff, wtLim);
00105    aveWt = MCeff*wtLim;
00106    cout<< "00000000000000000000000000000000000000000000000000000000000000000000000"<<endl;
00107    cout<< "00 -->wtLim: No_evt ="<<fNent<<"   <Wt> = "<<aveWt<<"  wtLim=  "<<wtLim<<endl;
00108    cout<< "00 -->wtLim: For eps = "<<eps  <<"    EFFICIENCY <Wt>/wtLim= "<<MCeff<<endl;
00109    cout<< "00000000000000000000000000000000000000000000000000000000000000000000000"<<endl;
00110 }
00111 
00112 //_________________________________________________________________________________
00113 void TFoamMaxwt::GetMCeff(Double_t eps, Double_t &MCeff, Double_t &wtLim)
00114 {
00115 // Calculates Efficiency= aveWt/wtLim for a given tolerance level epsilon<<1
00116 // using information stored in two histograms.
00117 // To be called at the end of the MC run.
00118 
00119    Int_t ib,ibX;
00120    Double_t lowEdge,bin,bin1;
00121    Double_t aveWt, aveWt1;
00122 
00123    fWtHst1->Print();
00124    fWtHst2->Print();
00125 
00126 // Convention on bin-numbering: nb=1 for 1-st bin, underflow nb=0, overflow nb=Nb+1
00127    Double_t sum   = 0.0;
00128    Double_t sumWt = 0.0;
00129    for(ib=0;ib<=fnBin+1;ib++) {
00130       sum   += fWtHst1->GetBinContent(ib);
00131       sumWt += fWtHst2->GetBinContent(ib);
00132    }
00133    if( (sum == 0.0) || (sumWt == 0.0) ) {
00134       cout<<"TFoamMaxwt::Make: zero content of histogram !!!,sum,sumWt ="<<sum<<sumWt<<endl;
00135    }
00136    aveWt = sumWt/sum;
00137    //--------------------------------------
00138    for( ibX=fnBin+1; ibX>0; ibX--) {
00139       lowEdge = (ibX-1.0)*fwmax/fnBin;
00140       sum   = 0.0;
00141       sumWt = 0.0;
00142       for( ib=0; ib<=fnBin+1; ib++) {
00143          bin  = fWtHst1->GetBinContent(ib);
00144          bin1 = fWtHst2->GetBinContent(ib);
00145          if(ib >= ibX) bin1=lowEdge*bin;
00146          sum   += bin;
00147          sumWt += bin1;
00148       }
00149       aveWt1 = sumWt/sum;
00150       if( TMath::Abs(1.0-aveWt1/aveWt) > eps ) break;
00151    }
00152    //---------------------------
00153    if(ibX == (fnBin+1) ) {
00154       wtLim = 1.0e200;
00155       MCeff   = 0.0;
00156       cout<< "+++++ wtLim undefined. Higher uper limit in histogram"<<endl;
00157    } else if( ibX == 1) {
00158       wtLim = 0.0;
00159       MCeff   =-1.0;
00160       cout<< "+++++ wtLim undefined. Lower uper limit or more bins "<<endl;
00161    } else {
00162       wtLim= (ibX)*fwmax/fnBin; // We over-estimate wtLim, under-estimate MCeff
00163       MCeff  = aveWt/wtLim;
00164    }
00165 }
00166 ///////////////////////////////////////////////////////////////////////////////
00167 //                                                                           //
00168 //      End of    Class  TFoamMaxwt                                          //
00169 //                                                                           //
00170 ///////////////////////////////////////////////////////////////////////////////

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