BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtHis2F Class Reference

#include <EvtHis2F.hh>

Inheritance diagram for EvtHis2F:

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)

Detailed Description

Definition at line 29 of file EvtHis2F.hh.

Constructor & Destructor Documentation

◆ EvtHis2F()

EvtHis2F::EvtHis2F ( )
inline

Definition at line 32 of file EvtHis2F.hh.

32{}

◆ ~EvtHis2F()

EvtHis2F::~EvtHis2F ( )
virtual

Definition at line 24 of file EvtHis2F.cc.

24{}

Member Function Documentation

◆ AR() [1/2]

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

Definition at line 199 of file EvtHis2F.cc.

199 { // accept: true. reject: false
200 // call this function after initiation and filling MC histogram HMC
201 bool accept;
202 accept = false;
203 if ( xmass2 < xlow || xmass2 > xup || ymass2 < ylow || ymass2 > yup ) { return accept; }
204 double zvalue = getZvalue( xmass2, ymass2 );
205 double ratio = zvalue / zmax;
206 double rndm = EvtRandom::Flat( 0.0, 1.0 );
207 // std::cout<<"zvalue, zmax , rndm "<<zvalue<<" "<<zmax<<" "<<rndm<<std::endl;
208 if ( ratio > rndm ) { accept = true; }
209 else { accept = false; }
210 return accept;
211}
double getZvalue(double m12_square, double m13_square)
Definition EvtHis2F.cc:190
static double Flat()
Definition EvtRandom.cc:69

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

◆ AR() [2/2]

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

Definition at line 213 of file EvtHis2F.cc.

214 { // accept: true. reject: false
215 // call this function after initiation and filling MC histogram HMC
216 bool accept;
217 accept = false;
218 setBins( h2 );
219 if ( xmass2 < xlow || xmass2 > xup || ymass2 < ylow || ymass2 > yup ) { return accept; }
220 int xbin = h2->GetXaxis()->FindBin( xmass2 );
221 int ybin = h2->GetYaxis()->FindBin( ymass2 );
222 int xybin = h2->GetBin( xbin, ybin );
223 double zvalue = h2->GetBinContent( xybin );
224
225 double ratio = zvalue / zmax;
226 double rndm = EvtRandom::Flat( 0.0, 1.0 );
227 // std::cout<<"zvalue, zmax , rndm "<<zvalue<<" "<<zmax<<" "<<rndm<<std::endl;
228 if ( ratio > rndm ) { accept = true; }
229 else { accept = false; }
230 return accept;
231}
void setBins(TH2F *h2)
Definition EvtHis2F.cc:144

◆ getBINSx()

int EvtHis2F::getBINSx ( )
inline

Definition at line 64 of file EvtHis2F.hh.

64{ return BINSx; }

◆ getBINSy()

int EvtHis2F::getBINSy ( )
inline

Definition at line 65 of file EvtHis2F.hh.

65{ return BINSy; }

◆ getFile()

const char * EvtHis2F::getFile ( )

Definition at line 26 of file EvtHis2F.cc.

26{ return datafile; }

Referenced by init().

◆ getHdt()

TH2F * EvtHis2F::getHdt ( )

Definition at line 52 of file EvtHis2F.cc.

52{ return HDATA; }

◆ getHmc()

TH2F * EvtHis2F::getHmc ( )

Definition at line 50 of file EvtHis2F.cc.

50{ return HMC; }

◆ getHTitle()

const char * EvtHis2F::getHTitle ( )

Definition at line 28 of file EvtHis2F.cc.

28{ return datatitle; }

Referenced by init().

◆ getHwt()

TH2F * EvtHis2F::getHwt ( )

Definition at line 54 of file EvtHis2F.cc.

54{ return HWT; }

Referenced by EvtNT3::init().

◆ getXlow()

double EvtHis2F::getXlow ( )
inline

Definition at line 66 of file EvtHis2F.hh.

66{ return xlow; }

◆ getXup()

double EvtHis2F::getXup ( )
inline

Definition at line 68 of file EvtHis2F.hh.

68{ return xup; }

◆ getYlow()

double EvtHis2F::getYlow ( )
inline

Definition at line 67 of file EvtHis2F.hh.

67{ return ylow; }

◆ getYup()

double EvtHis2F::getYup ( )
inline

Definition at line 69 of file EvtHis2F.hh.

69{ return yup; }

◆ getZmax()

double EvtHis2F::getZmax ( )

Definition at line 188 of file EvtHis2F.cc.

188{ return zmax; }

Referenced by EvtNT3::init().

◆ getZvalue()

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

Definition at line 190 of file EvtHis2F.cc.

190 {
191 int xbin = HWT->GetXaxis()->FindBin( xmass2 );
192 int ybin = HWT->GetYaxis()->FindBin( ymass2 );
193 int xybin = HWT->GetBin( xbin, ybin );
194 double zvalue = HWT->GetBinContent( xybin );
195 // std::cout<<"xmass, ymass, zvalue "<<xmass2<<" "<<ymass2<<" "<<zvalue<<std::endl;
196 return zvalue;
197}

Referenced by AR().

◆ HFill()

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

Definition at line 121 of file EvtHis2F.cc.

121 {
122 HMC->Fill( m12s, m13s, 1.0 );
123 icount++;
124 // if(icount==100000)show(HMC);
125 if ( icount == 100000 ) showFrame( HMC );
126 return;
127}
void showFrame(TH2F *h2)
Definition EvtHis2F.cc:74

◆ HReweight()

void EvtHis2F::HReweight ( )

Definition at line 129 of file EvtHis2F.cc.

129 { // reweight Dalitz Plot after get HDATA and HMC
130 std::cout << "Now I reweight the MC to the data Dalitz Plot" << std::endl;
131 double ndata = HDATA->Integral();
132 std::cout << "Number events in Dalitz plot = " << ndata << std::endl;
133 double nmc = HMC->Integral();
134 std::cout << "Number of events in HMC reweight= " << nmc << std::endl;
135 HWT->Divide( HDATA, HMC );
136 ndata = HWT->Integral();
137 std::cout << "Number of events after reweight= " << ndata << std::endl;
138 double norm = nmc / ndata;
139 HWT->Scale( norm );
140 nmc = HMC->Integral();
141 std::cout << "Number of events in HMC after scale= " << nmc << std::endl;
142}

Referenced by init().

◆ init() [1/2]

void EvtHis2F::init ( )

Definition at line 89 of file EvtHis2F.cc.

89 {
90 const char* dtfile = getFile(); // get data root file as input
91 const char* htitle = getHTitle(); // get Dalitz plot title
92
93 TFile dataf( dtfile );
94 HDATA = (TH2F*)dataf.Get( htitle );
95
96 double ndata = HDATA->Integral();
97 std::cout << "Number events in HDATA= " << ndata << std::endl;
98
99 TAxis* Xaxis = HDATA->GetXaxis();
100 TAxis* Yaxis = HDATA->GetYaxis();
101
102 BINSx = Xaxis->GetLast();
103 BINSy = Yaxis->GetLast();
104
105 xlow = Xaxis->GetBinLowEdge( 1 );
106 ylow = Yaxis->GetBinLowEdge( 1 );
107 xup = Xaxis->GetBinUpEdge( BINSx );
108 yup = Yaxis->GetBinUpEdge( BINSy );
109
110 std::cout << "BINSx,xlow,xup,BINSy,ylow,yup: " << BINSx << " ," << xlow << ", " << xup
111 << ", " << BINSy << ", " << ylow << ", " << yup << std::endl;
112 HMC = new TH2F( "myHMC", "", BINSx, xlow, xup, BINSy, ylow, yup );
113 HWT = new TH2F( "myHWT", "", BINSx, xlow, xup, BINSy, ylow, yup );
114 HMC->SetDirectory( 0 );
115 HWT->SetDirectory( 0 );
116 HDATA->SetDirectory( 0 );
117 icount = 0;
118 return;
119}
const char * getHTitle()
Definition EvtHis2F.cc:28
const char * getFile()
Definition EvtHis2F.cc:26

Referenced by EvtNT3::init().

◆ init() [2/2]

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

Definition at line 233 of file EvtHis2F.cc.

233 {
234 setHmc( hmc );
235 setHdt( hdt );
236 setHwt( hwt );
237 HReweight();
238 setZmax();
239}
void HReweight()
Definition EvtHis2F.cc:129
void setHwt(TH2F *H2)
Definition EvtHis2F.hh:50
void setHdt(TH2F *H2)
Definition EvtHis2F.cc:45
void setHmc(TH2F *H2)
Definition EvtHis2F.cc:40
void setZmax()
Definition EvtHis2F.cc:171

◆ setBins()

void EvtHis2F::setBins ( TH2F * h2)

Definition at line 144 of file EvtHis2F.cc.

144 {
145 TAxis* Xaxis = h2->GetXaxis();
146 TAxis* Yaxis = h2->GetYaxis();
147 BINSx = Xaxis->GetLast();
148 BINSy = Yaxis->GetLast();
149 xlow = Xaxis->GetBinLowEdge( 1 );
150 ylow = Yaxis->GetBinLowEdge( 1 );
151 xup = Xaxis->GetBinUpEdge( BINSx );
152 yup = Yaxis->GetBinUpEdge( BINSy );
153}

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

◆ setBINSx()

void EvtHis2F::setBINSx ( int bx)
inline

Definition at line 56 of file EvtHis2F.hh.

56{ BINSx = bx; }

◆ setBINSy()

void EvtHis2F::setBINSy ( int by)
inline

Definition at line 57 of file EvtHis2F.hh.

57{ BINSy = by; }

◆ setFile()

void EvtHis2F::setFile ( const char * dtfile)

Definition at line 30 of file EvtHis2F.cc.

30 {
31 datafile = dtfile;
32 return;
33}

◆ setHdt()

void EvtHis2F::setHdt ( TH2F * H2)

Definition at line 45 of file EvtHis2F.cc.

45 {
46 HDATA = (TH2F*)H2->Clone( "HDATA" );
47 return;
48}

Referenced by init().

◆ setHmc()

void EvtHis2F::setHmc ( TH2F * H2)

Definition at line 40 of file EvtHis2F.cc.

40 {
41 HMC = (TH2F*)H2->Clone( "HMC" );
42 return;
43}

Referenced by init().

◆ setHTitle()

void EvtHis2F::setHTitle ( const char * htitle)

Definition at line 35 of file EvtHis2F.cc.

35 {
36 datatitle = htitle;
37 return;
38}

◆ setHwt()

void EvtHis2F::setHwt ( TH2F * H2)
inline

Definition at line 50 of file EvtHis2F.hh.

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

Referenced by init().

◆ setXlow()

void EvtHis2F::setXlow ( double xl)
inline

Definition at line 58 of file EvtHis2F.hh.

58{ xlow = xl; }

◆ setXup()

void EvtHis2F::setXup ( double xu)
inline

Definition at line 59 of file EvtHis2F.hh.

59{ xup = xu; }

◆ setYlow()

void EvtHis2F::setYlow ( double yl)
inline

Definition at line 60 of file EvtHis2F.hh.

60{ ylow = yl; }

◆ setYup()

void EvtHis2F::setYup ( double yu)
inline

Definition at line 61 of file EvtHis2F.hh.

61{ yup = yu; }

◆ setZmax() [1/2]

void EvtHis2F::setZmax ( )

Definition at line 171 of file EvtHis2F.cc.

171 {
172 double ymax = 0;
173 setBins( HWT );
174
175 for ( int dx = 1; dx <= BINSx; dx++ )
176 {
177 for ( int dy = 1; dy <= BINSy; dy++ )
178 {
179 int nxy = HWT->GetBin( dx, dy );
180 double ncell = HWT->GetBinContent( nxy );
181 if ( ncell <= 0 ) continue;
182 if ( ncell > ymax ) { ymax = ncell; }
183 }
184 }
185 zmax = ymax;
186}

Referenced by init().

◆ setZmax() [2/2]

void EvtHis2F::setZmax ( TH2F * H2)

Definition at line 154 of file EvtHis2F.cc.

154 {
155 double ymax = 0;
156 setBins( hwt );
157
158 for ( int dx = 1; dx <= BINSx; dx++ )
159 {
160 for ( int dy = 1; dy <= BINSy; dy++ )
161 {
162 int nxy = hwt->GetBin( dx, dy );
163 double ncell = hwt->GetBinContent( nxy );
164 if ( ncell <= 0 ) continue;
165 if ( ncell > ymax ) { ymax = ncell; }
166 }
167 }
168 zmax = ymax;
169}

◆ show()

void EvtHis2F::show ( TH2F * h2)

Definition at line 56 of file EvtHis2F.cc.

56 {
57 TAxis* Xaxis = h2->GetXaxis();
58 TAxis* Yaxis = h2->GetYaxis();
59
60 int bx = Xaxis->GetLast();
61 int by = Yaxis->GetLast();
62
63 for ( int i = 1; i < bx + 1; i++ )
64 {
65 for ( int j = 1; j < by + 1; j++ )
66 {
67 int ij = h2->GetBin( i, j );
68 double xij = h2->GetBinContent( ij );
69 std::cout << "i,j,cell,value " << i << " " << j << " " << ij << " " << xij << std::endl;
70 }
71 }
72}

◆ showFrame()

void EvtHis2F::showFrame ( TH2F * h2)

Definition at line 74 of file EvtHis2F.cc.

74 {
75 TAxis* Xaxis = h2->GetXaxis();
76 TAxis* Yaxis = h2->GetYaxis();
77
78 int bx = Xaxis->GetLast();
79 int by = Yaxis->GetLast();
80
81 double x1 = Xaxis->GetBinLowEdge( 1 );
82 double y1 = Yaxis->GetBinLowEdge( 1 );
83 double x2 = Xaxis->GetBinUpEdge( bx );
84 double y2 = Yaxis->GetBinUpEdge( by );
85 std::cout << "N_x, x_low, x_up; N_y, y_low, y_up " << bx << " " << x1 << " " << x2 << " "
86 << by << " " << y1 << " " << y2 << std::endl;
87}

Referenced by HFill().


The documentation for this class was generated from the following files: