EvtHis2F Class Reference

#include <EvtHis2F.hh>

Inheritance diagram for EvtHis2F:

EvtNT3 List of all members.

Public Member Functions

 EvtHis2F ()
virtual ~EvtHis2F ()
void init ()
void setFile (const char *dtfile)
void setHTitle (const char *htitle)
const char * getFile ()
const char * getHTitle ()
void HFill (double xmass2, double ymass2)
double getZvalue (double m12_square, double m13_square)
void HReweight ()
bool AR (double xmass2, double ymass2)
bool AR (double xmass2, double ymass2, double zmax, TH2F *h2)
void setZmax ()
void setZmax (TH2F *H2)
double getZmax ()
void setHmc (TH2F *H2)
void setHdt (TH2F *H2)
void setHwt (TH2F *H2)
TH2F * getHmc ()
TH2F * getHdt ()
TH2F * getHwt ()
void setBINSx (int bx)
void setBINSy (int by)
void setXlow (double xl)
void setXup (double xu)
void setYlow (double yl)
void setYup (double yu)
void setBins (TH2F *h2)
int getBINSx ()
int getBINSy ()
double getXlow ()
double getYlow ()
double getXup ()
double getYup ()
void show (TH2F *h2)
void showFrame (TH2F *h2)
void init (TH2F *hmc, TH2F *hdt, TH2F *hwt)

Private Attributes

TH2F * HDATA
TH2F * HMC
TH2F * HWT
TFile * dataf
double zmax
int BINSx
int BINSy
double xlow
double xup
double ylow
double yup
const char * datafile
const char * datatitle
TDirectory dir
int icount
double max1
double max2
double max3
TH2F * WT1
TH2F * WT2
TH2F * WT3

Detailed Description

Definition at line 29 of file EvtHis2F.hh.


Constructor & Destructor Documentation

EvtHis2F::EvtHis2F (  )  [inline]

Definition at line 33 of file EvtHis2F.hh.

00033 {}

EvtHis2F::~EvtHis2F (  )  [virtual]

Definition at line 25 of file EvtHis2F.cc.

00025 {}


Member Function Documentation

bool EvtHis2F::AR ( double  xmass2,
double  ymass2,
double  zmax,
TH2F *  h2 
)

Definition at line 218 of file EvtHis2F.cc.

References EvtRandom::Flat(), setBins(), xup, and yup.

00218                                                                     {//accept: true. reject: false
00219   //call this function after initiation and filling MC histogram HMC
00220   bool accept;
00221   accept=false;
00222   setBins(h2);
00223   if(xmass2<xlow || xmass2>xup || ymass2<ylow || ymass2>yup) {return accept;}
00224   int xbin = h2->GetXaxis()->FindBin(xmass2);
00225   int ybin = h2->GetYaxis()->FindBin(ymass2);
00226   int xybin= h2->GetBin(xbin,ybin);
00227   double zvalue=h2->GetBinContent(xybin);
00228 
00229   double ratio=zvalue/zmax;
00230   double rndm=EvtRandom::Flat(0.0, 1.0);  
00231   //  std::cout<<"zvalue, zmax , rndm "<<zvalue<<" "<<zmax<<" "<<rndm<<std::endl;
00232   if(ratio > rndm){accept=true;} else {accept=false;}
00233   return accept;
00234 }

bool EvtHis2F::AR ( double  xmass2,
double  ymass2 
)

Definition at line 205 of file EvtHis2F.cc.

References EvtRandom::Flat(), getZvalue(), xup, yup, and zmax.

Referenced by EvtNT3::AR1(), EvtNT3::AR2(), EvtNT3::AR3(), and EvtDecay::SuperBody3decay_judge().

00205                                              {//accept: true. reject: false
00206   //call this function after initiation and filling MC histogram HMC
00207   bool accept;
00208   accept=false;
00209   if(xmass2<xlow || xmass2>xup || ymass2<ylow || ymass2>yup) {return accept;}
00210   double zvalue=getZvalue(xmass2,ymass2);
00211   double ratio=zvalue/zmax;
00212   double rndm=EvtRandom::Flat(0.0, 1.0);  
00213   //  std::cout<<"zvalue, zmax , rndm "<<zvalue<<" "<<zmax<<" "<<rndm<<std::endl;
00214   if(ratio > rndm){accept=true;} else {accept=false;}
00215   return accept;
00216 }

int EvtHis2F::getBINSx (  )  [inline]

Definition at line 64 of file EvtHis2F.hh.

References BINSx.

00064 {return BINSx;}

int EvtHis2F::getBINSy (  )  [inline]

Definition at line 65 of file EvtHis2F.hh.

References BINSy.

00065 {return BINSy;}

const char * EvtHis2F::getFile (  ) 

Definition at line 27 of file EvtHis2F.cc.

References datafile.

Referenced by init().

00027                              {
00028   return datafile;
00029 }

TH2F * EvtHis2F::getHdt (  ) 

Definition at line 59 of file EvtHis2F.cc.

References HDATA.

00059                        {
00060   return HDATA;
00061 }

TH2F * EvtHis2F::getHmc (  ) 

Definition at line 55 of file EvtHis2F.cc.

References HMC.

00055                        {
00056   return HMC;
00057 }

const char * EvtHis2F::getHTitle (  ) 

Definition at line 31 of file EvtHis2F.cc.

References datatitle.

Referenced by init().

00031                                {
00032   return datatitle;
00033 }

TH2F * EvtHis2F::getHwt (  ) 

Definition at line 63 of file EvtHis2F.cc.

References HWT.

Referenced by EvtNT3::init().

00063                        {
00064   return HWT;
00065 }

double EvtHis2F::getXlow (  )  [inline]

Definition at line 66 of file EvtHis2F.hh.

References xlow.

Referenced by EvtDecay::SuperBody3decay_make().

00066 {return xlow;}

double EvtHis2F::getXup (  )  [inline]

Definition at line 68 of file EvtHis2F.hh.

References xup.

Referenced by EvtDecay::SuperBody3decay_make().

00068 {return xup;}

double EvtHis2F::getYlow (  )  [inline]

Definition at line 67 of file EvtHis2F.hh.

References ylow.

Referenced by EvtDecay::SuperBody3decay_make().

00067 {return ylow;}

double EvtHis2F::getYup (  )  [inline]

Definition at line 69 of file EvtHis2F.hh.

References yup.

Referenced by EvtDecay::SuperBody3decay_make().

00069 {return yup;}

double EvtHis2F::getZmax (  ) 

Definition at line 192 of file EvtHis2F.cc.

References zmax.

Referenced by EvtNT3::init().

00192                          {
00193   return zmax;
00194 }

double EvtHis2F::getZvalue ( double  m12_square,
double  m13_square 
)

Definition at line 196 of file EvtHis2F.cc.

References HWT.

Referenced by AR().

00196                                                       {
00197   int xbin = HWT->GetXaxis()->FindBin(xmass2);
00198   int ybin = HWT->GetYaxis()->FindBin(ymass2);
00199   int xybin= HWT->GetBin(xbin,ybin);
00200   double zvalue=HWT->GetBinContent(xybin);
00201 //  std::cout<<"xmass, ymass, zvalue "<<xmass2<<" "<<ymass2<<" "<<zvalue<<std::endl;
00202   return zvalue;
00203 }

void EvtHis2F::HFill ( double  xmass2,
double  ymass2 
)

Definition at line 129 of file EvtHis2F.cc.

References HMC, icount, and showFrame().

Referenced by EvtDecay::SuperBody3decay_make().

00129                                             {
00130   HMC->Fill(m12s,m13s,1.0);
00131   icount++;
00132   //  if(icount==100000)show(HMC);
00133   if(icount==100000) showFrame(HMC);
00134   return;
00135 }

void EvtHis2F::HReweight (  ) 

Definition at line 137 of file EvtHis2F.cc.

References HDATA, HMC, and HWT.

Referenced by init(), and EvtDecay::SuperBody3decay_make().

00137                          {  //reweight Dalitz Plot after get HDATA and HMC
00138   std::cout<<"Now I reweight the MC to the data Dalitz Plot"<<std::endl;
00139   double ndata=HDATA->Integral();
00140   std::cout<<"Number events in Dalitz plot = "<<ndata<<std::endl;
00141   double nmc=HMC->Integral();
00142   std::cout<<"Number of events in HMC reweight= "<<nmc<<std::endl;
00143   HWT->Divide(HDATA,HMC);
00144   ndata=HWT->Integral();
00145   std::cout<<"Number of events after reweight= "<<ndata<<std::endl;
00146   double norm=nmc/ndata;
00147   HWT->Scale(norm);
00148   nmc=HMC->Integral();
00149   std::cout<<"Number of events in HMC after scale= "<<nmc<<std::endl;
00150 }

void EvtHis2F::init ( TH2F *  hmc,
TH2F *  hdt,
TH2F *  hwt 
)

Definition at line 236 of file EvtHis2F.cc.

References HReweight(), setHdt(), setHmc(), setHwt(), and setZmax().

00236                                                   {
00237   setHmc(hmc);
00238   setHdt(hdt);
00239   setHwt(hwt);
00240   HReweight();
00241   setZmax();
00242 }

void EvtHis2F::init (  ) 

Reimplemented in EvtNT3.

Definition at line 98 of file EvtHis2F.cc.

References BINSx, BINSy, dataf, getFile(), getHTitle(), HDATA, HMC, HWT, icount, xlow, xup, ylow, and yup.

Referenced by EvtNT3::init(), and EvtDecay::initialize().

00098                     {
00099   const char* dtfile =getFile();  // get data root file as input
00100   const char* htitle =getHTitle();  // get Dalitz plot title
00101 
00102   TFile dataf(dtfile);
00103   HDATA = (TH2F*)dataf.Get(htitle);
00104 
00105   double ndata=HDATA->Integral();
00106   std::cout<<"Number events in HDATA= "<<ndata<<std::endl;
00107 
00108   TAxis* Xaxis=HDATA->GetXaxis();
00109   TAxis* Yaxis=HDATA->GetYaxis();
00110  
00111   BINSx  =Xaxis->GetLast();  
00112   BINSy  =Yaxis->GetLast();
00113   
00114   xlow=Xaxis->GetBinLowEdge(1);
00115   ylow=Yaxis->GetBinLowEdge(1);
00116   xup =Xaxis->GetBinUpEdge(BINSx);
00117   yup =Yaxis->GetBinUpEdge(BINSy);
00118 
00119   std::cout<<"BINSx,xlow,xup,BINSy,ylow,yup: "<<BINSx<<" ,"<<xlow<<", "<<xup<<", "<<BINSy<<", "<<ylow<<", "<<yup<<std::endl;
00120   HMC = new TH2F("myHMC","",BINSx,xlow,xup,BINSy,ylow,yup);
00121   HWT = new TH2F("myHWT","",BINSx,xlow,xup,BINSy,ylow,yup);
00122   HMC ->SetDirectory(0);
00123   HWT ->SetDirectory(0);
00124   HDATA->SetDirectory(0);
00125   icount=0;
00126   return;
00127 }

void EvtHis2F::setBins ( TH2F *  h2  ) 

Definition at line 152 of file EvtHis2F.cc.

References BINSx, BINSy, xlow, xup, ylow, and yup.

Referenced by AR(), and setZmax().

00152                                {
00153   TAxis* Xaxis=h2->GetXaxis();
00154   TAxis* Yaxis=h2->GetYaxis();
00155   BINSx  =Xaxis->GetLast();  
00156   BINSy  =Yaxis->GetLast();
00157   xlow=Xaxis->GetBinLowEdge(1);
00158   ylow=Yaxis->GetBinLowEdge(1);
00159   xup =Xaxis->GetBinUpEdge(BINSx);
00160   yup =Yaxis->GetBinUpEdge(BINSy);
00161 }

void EvtHis2F::setBINSx ( int  bx  )  [inline]

Definition at line 56 of file EvtHis2F.hh.

References BINSx.

00056 {BINSx=bx;}

void EvtHis2F::setBINSy ( int  by  )  [inline]

Definition at line 57 of file EvtHis2F.hh.

References BINSy.

00057 {BINSy=by;}

void EvtHis2F::setFile ( const char *  dtfile  ) 

Definition at line 35 of file EvtHis2F.cc.

References datafile.

Referenced by EvtDecay::initialize().

00035                                         {
00036   datafile=dtfile;
00037   return;
00038 }

void EvtHis2F::setHdt ( TH2F *  H2  ) 

Definition at line 50 of file EvtHis2F.cc.

References HDATA.

Referenced by init().

00050                                {
00051   HDATA=(TH2F*) H2->Clone("HDATA");
00052   return;
00053 }

void EvtHis2F::setHmc ( TH2F *  H2  ) 

Definition at line 45 of file EvtHis2F.cc.

References HMC.

Referenced by init().

00045                                {
00046   HMC=(TH2F*) H2->Clone("HMC");
00047   return;
00048 }

void EvtHis2F::setHTitle ( const char *  htitle  ) 

Definition at line 40 of file EvtHis2F.cc.

References datatitle.

Referenced by EvtDecay::initialize().

00040                                           {
00041   datatitle=htitle;
00042   return;
00043 }

void EvtHis2F::setHwt ( TH2F *  H2  )  [inline]

Definition at line 50 of file EvtHis2F.hh.

References HWT.

Referenced by init().

00050 { HWT=(TH2F*) H2->Clone("HWT");}

void EvtHis2F::setXlow ( double  xl  )  [inline]

Definition at line 58 of file EvtHis2F.hh.

References xlow.

00058 {xlow=xl;}

void EvtHis2F::setXup ( double  xu  )  [inline]

Definition at line 59 of file EvtHis2F.hh.

References xup.

00059 {xup=xu;}

void EvtHis2F::setYlow ( double  yl  )  [inline]

Definition at line 60 of file EvtHis2F.hh.

References ylow.

00060 {ylow=yl;}

void EvtHis2F::setYup ( double  yu  )  [inline]

Definition at line 61 of file EvtHis2F.hh.

References yup.

00061 {yup=yu;}

void EvtHis2F::setZmax ( TH2F *  H2  ) 

Definition at line 162 of file EvtHis2F.cc.

References BINSx, BINSy, setBins(), and zmax.

00162                                 {
00163   double ymax=0;
00164   setBins(hwt);
00165 
00166    for(int dx=1;dx <= BINSx;dx++){
00167      for(int dy=1;dy <= BINSy;dy++){
00168        int nxy=hwt->GetBin(dx,dy);
00169        double ncell=hwt->GetBinContent(nxy);
00170        if(ncell<=0) continue;
00171        if(ncell>ymax){ymax=ncell;}
00172      }
00173    }
00174    zmax=ymax;
00175 }

void EvtHis2F::setZmax (  ) 

Definition at line 177 of file EvtHis2F.cc.

References BINSx, BINSy, HWT, setBins(), and zmax.

Referenced by init(), and EvtDecay::SuperBody3decay_make().

00177                        {
00178   double ymax=0;
00179   setBins(HWT);
00180 
00181    for(int dx=1;dx <= BINSx;dx++){
00182      for(int dy=1;dy <= BINSy;dy++){
00183        int nxy=HWT->GetBin(dx,dy);
00184        double ncell=HWT->GetBinContent(nxy);
00185        if(ncell<=0) continue;
00186        if(ncell>ymax){ymax=ncell;}
00187      }
00188    }
00189    zmax=ymax;
00190 }

void EvtHis2F::show ( TH2F *  h2  ) 

Definition at line 67 of file EvtHis2F.cc.

References genRecEmupikp::i, and ganga-rec::j.

00067                             {
00068   TAxis* Xaxis=h2->GetXaxis();
00069   TAxis* Yaxis=h2->GetYaxis();
00070  
00071   int bx  =Xaxis->GetLast();  
00072   int by  =Yaxis->GetLast();
00073 
00074   for(int i=1;i<bx+1;i++){
00075     for(int j=1;j<by+1;j++){
00076       int ij=h2->GetBin(i,j);
00077       double xij=h2->GetBinContent(ij);
00078       std::cout<<"i,j,cell,value "<<i<<" "<<j<<" "<<ij<<" "<<xij<<std::endl;
00079     }
00080   }
00081 }

void EvtHis2F::showFrame ( TH2F *  h2  ) 

Definition at line 83 of file EvtHis2F.cc.

Referenced by HFill().

00083                                  {
00084   TAxis* Xaxis=h2->GetXaxis();
00085   TAxis* Yaxis=h2->GetYaxis();
00086  
00087   int bx  =Xaxis->GetLast();  
00088   int by  =Yaxis->GetLast();
00089   
00090   double x1 =Xaxis->GetBinLowEdge(1);
00091   double y1 =Yaxis->GetBinLowEdge(1);
00092   double x2 =Xaxis->GetBinUpEdge(bx);
00093   double y2 =Yaxis->GetBinUpEdge(by);
00094   std::cout<<"N_x, x_low, x_up; N_y, y_low, y_up "
00095            <<bx<<" "<<x1<<" "<<x2<<" "<<by<<" "<<y1<<" "<<y2<<std::endl;
00096 }


Member Data Documentation

int EvtHis2F::BINSx [private]

Definition at line 80 of file EvtHis2F.hh.

Referenced by getBINSx(), init(), setBins(), setBINSx(), and setZmax().

int EvtHis2F::BINSy [private]

Definition at line 80 of file EvtHis2F.hh.

Referenced by getBINSy(), init(), setBins(), setBINSy(), and setZmax().

TFile* EvtHis2F::dataf [private]

Reimplemented in EvtNT3.

Definition at line 78 of file EvtHis2F.hh.

Referenced by init().

const char* EvtHis2F::datafile [private]

Reimplemented in EvtNT3.

Definition at line 82 of file EvtHis2F.hh.

Referenced by getFile(), and setFile().

const char * EvtHis2F::datatitle [private]

Definition at line 82 of file EvtHis2F.hh.

Referenced by getHTitle(), and setHTitle().

TDirectory EvtHis2F::dir [private]

Definition at line 83 of file EvtHis2F.hh.

TH2F* EvtHis2F::HDATA [private]

Definition at line 77 of file EvtHis2F.hh.

Referenced by getHdt(), HReweight(), init(), and setHdt().

TH2F * EvtHis2F::HMC [private]

Definition at line 77 of file EvtHis2F.hh.

Referenced by getHmc(), HFill(), HReweight(), init(), and setHmc().

TH2F * EvtHis2F::HWT [private]

Definition at line 77 of file EvtHis2F.hh.

Referenced by getHwt(), getZvalue(), HReweight(), init(), setHwt(), and setZmax().

int EvtHis2F::icount [private]

Definition at line 84 of file EvtHis2F.hh.

Referenced by HFill(), and init().

double EvtHis2F::max1 [private]

Reimplemented in EvtNT3.

Definition at line 85 of file EvtHis2F.hh.

double EvtHis2F::max2 [private]

Reimplemented in EvtNT3.

Definition at line 85 of file EvtHis2F.hh.

double EvtHis2F::max3 [private]

Reimplemented in EvtNT3.

Definition at line 85 of file EvtHis2F.hh.

TH2F* EvtHis2F::WT1 [private]

Reimplemented in EvtNT3.

Definition at line 86 of file EvtHis2F.hh.

TH2F * EvtHis2F::WT2 [private]

Reimplemented in EvtNT3.

Definition at line 86 of file EvtHis2F.hh.

TH2F * EvtHis2F::WT3 [private]

Reimplemented in EvtNT3.

Definition at line 86 of file EvtHis2F.hh.

double EvtHis2F::xlow [private]

Definition at line 81 of file EvtHis2F.hh.

Referenced by getXlow(), init(), setBins(), and setXlow().

double EvtHis2F::xup [private]

Definition at line 81 of file EvtHis2F.hh.

Referenced by AR(), getXup(), init(), setBins(), and setXup().

double EvtHis2F::ylow [private]

Definition at line 81 of file EvtHis2F.hh.

Referenced by getYlow(), init(), setBins(), and setYlow().

double EvtHis2F::yup [private]

Definition at line 81 of file EvtHis2F.hh.

Referenced by AR(), getYup(), init(), setBins(), and setYup().

double EvtHis2F::zmax [private]

Definition at line 79 of file EvtHis2F.hh.

Referenced by AR(), getZmax(), and setZmax().


Generated on Tue Nov 29 23:19:01 2016 for BOSS_7.0.2 by  doxygen 1.4.7