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

#include <EvtNT3.hh>

Inheritance diagram for EvtNT3:

Public Member Functions

void init ()
int calculateBins (int entries)
bool AR1 (double costheta, double mass)
bool AR2 (double costheta, double mass)
bool AR3 (double costheta, double mass)
void setDTfile (const char *name)
void setMCfile (const char *name)
Public Member Functions inherited from EvtHis2F
 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 25 of file EvtNT3.hh.

Member Function Documentation

◆ AR1()

bool EvtNT3::AR1 ( double costheta,
double mass )

Definition at line 141 of file EvtNT3.cc.

141 {
142 bool accept = false;
143 accept = EvtHis2F::AR( costheta, mass, max1, WT1 );
144 //--- debugging
145 // std::cout<<"Max_mass= "<<getZmax()<<std::endl;
146
147 return accept;
148}
double mass
bool AR(double xmass2, double ymass2)
Definition EvtHis2F.cc:199

◆ AR2()

bool EvtNT3::AR2 ( double costheta,
double mass )

Definition at line 150 of file EvtNT3.cc.

150 {
151 bool accept = false;
152 accept = EvtHis2F::AR( costheta, mass, max2, WT2 );
153
154 return accept;
155}

◆ AR3()

bool EvtNT3::AR3 ( double costheta,
double mass )

Definition at line 157 of file EvtNT3.cc.

157 {
158 bool accept = false;
159 accept = EvtHis2F::AR( costheta, mass, max3, WT3 );
160
161 return accept;
162}

◆ calculateBins()

int EvtNT3::calculateBins ( int entries)

Definition at line 133 of file EvtNT3.cc.

133 {
134 int bins, ncell;
135 ncell = 30; // at lease to require each cell to have 30 events
136 bins = entries / ncell / Ncos;
137 if ( bins > 100 ) { bins = 100; }
138 return bins;
139}
int bins[20]

Referenced by init().

◆ init()

void EvtNT3::init ( )

Definition at line 25 of file EvtNT3.cc.

25 {
26 max = 0;
27 Ncos = 20;
28 chainMC = new TChain( "mc" );
29 chainDT = new TChain( "data" );
30
31 chainMC->SetDirectory( 0 );
32 chainDT->SetDirectory( 0 );
33
34 chainMC->Add( mcfile );
35 chainDT->Add( datafile );
36 //--MC
37 chainMC->SetBranchAddress( "costheta1", &costheta1 );
38 chainMC->SetBranchAddress( "costheta2", &costheta2 );
39 chainMC->SetBranchAddress( "costheta3", &costheta3 );
40 chainMC->SetBranchAddress( "m12", &m12 );
41 chainMC->SetBranchAddress( "m13", &m13 );
42 chainMC->SetBranchAddress( "m23", &m23 );
43 //--Data
44 chainDT->SetBranchAddress( "costheta1", &costheta1 );
45 chainDT->SetBranchAddress( "costheta2", &costheta2 );
46 chainDT->SetBranchAddress( "costheta3", &costheta3 );
47 chainDT->SetBranchAddress( "m12", &m12 );
48 chainDT->SetBranchAddress( "m13", &m13 );
49 chainDT->SetBranchAddress( "m23", &m23 );
50
51 entriesMC = (Int_t)chainMC->GetEntries();
52 entriesDT = (Int_t)chainDT->GetEntries();
53
54 m12_low = chainDT->GetMinimum( "m12" );
55 m13_low = chainDT->GetMinimum( "m13" );
56 m23_low = chainDT->GetMinimum( "m23" );
57
58 m12_up = chainDT->GetMaximum( "m12" );
59 m13_up = chainDT->GetMaximum( "m13" );
60 m23_up = chainDT->GetMaximum( "m23" );
61
62 Int_t ny = calculateBins( entriesDT ); // bins in mass axisis (Y-axisis)
63
64 MC1 = new TH2F( "myMC1", "", Ncos, -1.0, 1.0, ny, m23_low, m23_up ); // costheta1 vs. m23
65 MC2 = new TH2F( "myMC2", "", Ncos, -1.0, 1.0, ny, m13_low, m13_up ); // costheta2 vs. m13
66 MC3 = new TH2F( "myMC3", "", Ncos, -1.0, 1.0, ny, m12_low, m12_up ); // costheta3 vs. m12
67
68 DT1 = new TH2F( "myDT1", "", Ncos, -1.0, 1.0, ny, m23_low, m23_up ); // costheta1 vs. m23
69 DT2 = new TH2F( "myDT2", "", Ncos, -1.0, 1.0, ny, m13_low, m13_up ); // costheta2 vs. m13
70 DT3 = new TH2F( "myDT3", "", Ncos, -1.0, 1.0, ny, m12_low, m12_up ); // costheta3 vs. m12
71
72 WT1 = new TH2F( "myWT1", "", Ncos, -1.0, 1.0, ny, m23_low, m23_up ); // costheta1 vs. m23
73 WT2 = new TH2F( "myWT2", "", Ncos, -1.0, 1.0, ny, m13_low, m13_up ); // costheta2 vs. m13
74 WT3 = new TH2F( "myWT3", "", Ncos, -1.0, 1.0, ny, m12_low, m12_up ); // costheta3 vs. m12
75
76 MC1->SetDirectory( 0 );
77 MC2->SetDirectory( 0 );
78 MC3->SetDirectory( 0 );
79
80 DT1->SetDirectory( 0 );
81 DT2->SetDirectory( 0 );
82 DT3->SetDirectory( 0 );
83
84 WT1->SetDirectory( 0 );
85 WT2->SetDirectory( 0 );
86 WT3->SetDirectory( 0 );
87
88 // filling MC histogram
89 for ( Int_t j = 0; j < entriesMC; j++ )
90 {
91 chainMC->GetEntry( j );
92 MC1->Fill( costheta1, m23 );
93 MC2->Fill( costheta2, m13 );
94 MC3->Fill( costheta3, m12 );
95 }
96
97 // filling data histogram
98 for ( Int_t j = 0; j < entriesDT; j++ )
99 {
100 chainDT->GetEntry( j );
101 DT1->Fill( costheta1, m23 );
102 DT2->Fill( costheta2, m13 );
103 DT3->Fill( costheta3, m12 );
104 }
105 //--------debugging ---------
106 // showFrame(MC1);
107 // showFrame(MC2);
108 // showFrame(MC3);
109
110 EvtHis2F::init( MC1, DT1, WT1 );
111 WT1 = EvtHis2F::getHwt();
112 max1 = EvtHis2F::getZmax();
113
114 EvtHis2F::init( MC2, DT2, WT2 );
115 WT2 = EvtHis2F::getHwt();
116 max2 = EvtHis2F::getZmax();
117
118 EvtHis2F::init( MC3, DT3, WT3 );
119 WT3 = EvtHis2F::getHwt();
120 max3 = EvtHis2F::getZmax();
121
122 //------------
123
124 /*
125 show(WT1);
126 std::cout<<"================================================="<<std::endl;
127 show(WT2);
128 std::cout<<"================================================="<<std::endl;
129 show(WT3);
130 */
131}
void init()
Definition EvtHis2F.cc:89
double getZmax()
Definition EvtHis2F.cc:188
TH2F * getHwt()
Definition EvtHis2F.cc:54
int calculateBins(int entries)
Definition EvtNT3.cc:133

◆ setDTfile()

void EvtNT3::setDTfile ( const char * name)
inline

Definition at line 37 of file EvtNT3.hh.

37{ datafile = name; }

◆ setMCfile()

void EvtNT3::setMCfile ( const char * name)
inline

Definition at line 38 of file EvtNT3.hh.

38{ mcfile = name; }

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