BOSS 7.1.1
BESIII Offline Software System
Loading...
Searching...
No Matches
TofDataSet.cxx
Go to the documentation of this file.
2#include <cmath>
3
5 initial();
6 return;
7}
8
9
11
12
14 m_run = 0;
15 m_event = -1;
16 m_tofid = -1;
17 m_strip = -1;
18 m_qleft = -99.0;
19 m_qright = -99.0;
20 m_tleft = -99.0;
21 m_tright = -99.0;
22 m_zrhit = -999.0;
23 m_dt = -999.0;
24 m_texp = -99.0;
25 m_path = -99.0;
26 m_phi = -99.0;
27 m_theta = -99.0;
28 m_p = -9.0;
29 m_t0 = -99.0;
30 m_q0 = -99.0;
31 m_hitcase = -1;
32 return;
33}
34
35
37 initial();
38 if( hit ) {
39 m_run = hit->run();
40 m_event = hit->event();
41 m_tofid = hit->mod();
42 m_strip = int( hit->sinTheta() );
43 m_qleft = hit->adc1();
44 m_qright = hit->adc2();
45 m_tleft = hit->tdc1();
46 m_tright = hit->tdc2();
47 m_zrhit = hit->zHit();
48 m_dt = ( hit->tdc2() - hit->tdc1() )/2.0;
49 m_texp = hit->tpred();
50 m_path = hit->path();
51 m_phi = hit->deltaPhi();
52 m_theta = hit->sinTheta();
53 m_p = hit->p();
54 m_hitcase = hit->qual();
55 }
56 return;
57}
58
59
61 initial();
62 if( hit ) {
63 m_run = hit->run();
64 m_event = hit->event();
65 m_tofid = hit->mod();
66 m_qleft = hit->adc();
67 m_tleft = hit->tdc();
68 m_zrhit = hit->rHit();
69 m_texp = hit->tpred();
70 m_path = hit->path();
71 m_phi = hit->deltaPhi();
72 m_theta = hit->cosTheta();
73 m_p = hit->p();
74 m_hitcase = hit->qual();
75 }
76 return;
77}
78
79
81 initial();
82 if( &one ) {
83 m_run = one.run;
84 m_event = one.event;
85 m_tofid = one.tofid;
86 m_strip = one.strip;
87 m_qleft = one.qleft;
88 m_tleft = one.tleft;
89 if( ( one.hitcase >= 0 && one.hitcase <=2 ) || ( one.hitcase==5 || one.hitcase==6 ) ) {
90 m_qright = one.qright;
91 m_tright = one.tright;
92 }
93 m_zrhit = one.zrhit;
94 m_dt = one.dt;
95 m_texp = one.texp;
96 m_path = one.path;
97 m_phi = one.phi;
98 m_theta = one.theta;
99 m_p = one.p;
100 m_hitcase = one.hitcase;
101 }
102 return;
103}
104
105
106void TofDataSet::setData( TTree* t, unsigned int isBarrel ) {
107 if( t && t->GetEntries()>0 ) {
108 rootRecord item;
109 t->SetBranchAddress( "run", &item.run );
110 t->SetBranchAddress( "event", &item.event );
111 t->SetBranchAddress( "tofid", &item.tofid );
112 if( NULL != t->FindBranch( "strip" ) ) {
113 t->SetBranchAddress( "strip", &item.strip );
114 }
115 t->SetBranchAddress( "qleft", &item.qleft );
116 t->SetBranchAddress( "qright", &item.qright );
117 t->SetBranchAddress( "tleft", &item.tleft );
118 t->SetBranchAddress( "tright", &item.tright );
119 t->SetBranchAddress( "zrhit", &item.zrhit );
120 if( NULL != t->FindBranch( "dt" ) ) {
121 t->SetBranchAddress( "dt", &item.dt );
122 }
123 t->SetBranchAddress( "texp", &item.texp );
124 t->SetBranchAddress( "path", &item.path );
125 t->SetBranchAddress( "phi", &item.phi );
126 t->SetBranchAddress( "theta", &item.theta );
127 t->SetBranchAddress( "p", &item.p );
128 t->SetBranchAddress( "hitcase", &item.hitcase );
129
130 for( unsigned int i=0; i<t->GetEntries(); i++ ) {
131 t->GetEntry(i);
132 if( isBarrel==1 && item.hitcase>=0 && item.hitcase<=2 ) {
133 Record *r = new Record( item );
134 if( r->cutBarrel() ) {
135 unsigned int tofID = item.tofid;
136 barrelData[tofID]->push_back(r);
137 }
138 else {
139 delete r;
140 }
141 }
142 else if( isBarrel==0 && ( item.hitcase==3 || item.hitcase==4 ) ) {
143 Record *r = new Record( item );
144 if( r->cutEndcap() ) {
145 unsigned int tofID = item.tofid;
146 endcapData[tofID]->push_back(r);
147 }
148 else {
149 delete r;
150 }
151 }
152 else if( isBarrel==2 && item.hitcase>=5 && item.hitcase<=6 ) {
153 Record *r = new Record( item );
154 if( r->cutEtf() ) {
155 unsigned int tofID = item.tofid;
156 unsigned int strip = item.strip;
157 unsigned int id = tofID*12 + strip;
158 etfData[id]->push_back(r);
159 }
160 else {
161 delete r;
162 }
163 }
164 }
165 }
166 else {
167 std::cerr << "Error: a invalid tree or a blank tree, When converting a tree to TofDataSet,exit" << std::endl;
168 }
169
170 return;
171}
172
173
175 for( unsigned int i=0; i<NBarrel; i++ ) {
176 barrelData[i] = new RecordSet;
177 }
178 for( unsigned int i=0; i<NEndcap; i++ ) {
179 endcapData[i] = new RecordSet;
180 }
181 for( unsigned int i=0; i<NEtf*NStrip; i++ ) {
182 etfData[i] = new RecordSet;
183 }
184 return;
185}
186
187
189 for( unsigned int i=0; i<NBarrel; i++ ) {
190 barrelData[i]->clear();
191 delete barrelData[i];
192 }
193 for( unsigned int i=0; i<NEndcap; i++ ) {
194 endcapData[i]->clear();
195 delete endcapData[i];
196 }
197 for( unsigned int i=0; i<NEtf*NStrip; i++ ) {
198 etfData[i]->clear();
199 delete etfData[i];
200 }
201 return;
202}
203
204
205void TofDataSet::setBarrelDataFiles( std::vector<std::string>& barrelFiles ) {
206 TChain* data_barrel = new TChain("btrk");
207 if( !data_barrel ) {
208 std::cerr << " tofcalgsec Error Msg: creating a tree[barrel] fails in TofDataSet()"<<std::endl;
209 throw "Error Msg: creating a tree fails in TofDataSet() ";
210 }
211 std::cout<<"begin reading barrel data file ... "<<std::endl;
212 try{
213 for( std::vector<std::string>::iterator it=barrelFiles.begin(); it!=barrelFiles.end(); it++ ) {
214 std::cout << " Add file : " << (*it) << std::endl;
215 data_barrel->Add( (*it).c_str() );
216 }
217 }
218 catch(...){
219 std::cerr << "tofcalgsec Error Msg : in TofDataSet::setDataFiles(std::vector<std::string>&) " << std::endl;
220 return;
221 }
222 setData( data_barrel, 1 );
223 delete data_barrel;
224 return;
225}
226
227
228void TofDataSet::setEndcapDataFiles( std::vector<std::string>& endcapFiles ) {
229 TChain* data_endcap = new TChain("etrk");
230 if( !data_endcap ) {
231 std::cerr << " tofcalgsec Error Msg: creating a tree[endcap] fails in TofDataSet()"<<std::endl;
232 throw "Error Msg: creating a tree fails in TofDataSet() ";
233 }
234 std::cout<<"begin reading endcap data file ... "<<std::endl;
235 try{
236 for( std::vector<std::string>::iterator it=endcapFiles.begin(); it!=endcapFiles.end(); it++ ) {
237 std::cout << " Add file : " << (*it) << std::endl;
238 data_endcap->Add( (*it).c_str() );
239 }
240 }
241 catch(...){
242 std::cerr << "tofcalgsec Error Msg : in TofDataSet::setDataFiles(std::vector<std::string>&) " << std::endl;
243 return;
244 }
245 setData( data_endcap, 0 );
246 delete data_endcap;
247 return;
248}
249
250
251void TofDataSet::setEtfDataFiles( std::vector<std::string>& etfFiles ) {
252 TChain* data_etf = new TChain("etf");
253 if( !data_etf ) {
254 std::cerr << " tofcalgsec Error Msg: creating a tree[etf] fails in TofDataSet()"<<std::endl;
255 throw "Error Msg: creating a tree fails in TofDataSet() ";
256 }
257 std::cout<<"begin reading etf data file ... "<<std::endl;
258 try{
259 for( std::vector<std::string>::iterator it=etfFiles.begin(); it!=etfFiles.end(); it++ ) {
260 std::cout << " Add file : " << (*it) << std::endl;
261 data_etf->Add( (*it).c_str() );
262 }
263 }
264 catch(...){
265 std::cerr << "tofcalgsec Error Msg : in TofDataSet::setDataFiles(std::vector<std::string>&) " << std::endl;
266 return;
267 }
268 setData( data_etf, 2 );
269 delete data_etf;
270 return;
271}
272
273
275 RecBTofCalHitCol::iterator iter = bhitcol.begin();
276 for( ; iter!=bhitcol.end(); iter++ ) {
277 int tofid = (*iter)->mod();
278 if( tofid<0 || tofid>175 ) continue;
279 if( fabs( (*iter)->dzHit() - 1.0 )>1.0e-6 ) continue;
280
281 Record *r = new Record( (*iter) );
282 if( r->cutBarrel() ) {
283 barrelData[tofid]->push_back(r);
284 }
285 else {
286 delete r;
287 }
288 }
289 return;
290}
291
292
294 RecETofCalHitCol::iterator iter = ehitcol.begin();
295 for( ; iter!=ehitcol.end(); iter++ ) {
296 int tofid = (*iter)->mod();
297 if( tofid<0 || tofid>95 ) continue;
298
299 Record *r = new Record( (*iter) );
300 if( r->cutEndcap() ) {
301 endcapData[tofid]->push_back(r);
302 }
303 else {
304 delete r;
305 }
306 }
307 return;
308}
309
310
312 RecBTofCalHitCol::iterator iter = bhitcol.begin();
313 for( ; iter!=bhitcol.end(); iter++ ) {
314 int tofid = (*iter)->mod();
315 int strip = int( (*iter)->sinTheta() );
316 if( tofid<0 || tofid>71 ) continue;
317 if( strip<0 || strip>11 ) continue;
318 if( fabs( (*iter)->dzHit() )>1.0e-6 ) continue;
319 unsigned int id = tofid*12 + strip;
320 Record *r = new Record( (*iter) );
321 if( r->cutEtf() ) {
322 etfData[id]->push_back(r);
323 }
324 else {
325 delete r;
326 }
327 }
328 return;
329}
330
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
ObjectVector< RecBTofCalHit > RecBTofCalHitCol
ObjectVector< RecETofCalHit > RecETofCalHitCol
std::vector< Record * > RecordSet
Definition TofDataSet.h:98
const unsigned int NStrip
Definition TofDataSet.h:15
const unsigned int NEtf
Definition TofDataSet.h:14
const unsigned int NBarrel
Definition TofDataSet.h:12
const unsigned int NEndcap
Definition TofDataSet.h:13
#define NULL
TTree * t
Definition binning.cxx:23
double tpred() const
double p() const
double zHit() const
double tdc1() const
int event() const
int mod() const
double adc1() const
double deltaPhi() const
double tdc2() const
int qual() const
double adc2() const
int run() const
double sinTheta() const
double path() const
double p() const
double tpred() const
double deltaPhi() const
int run() const
double path() const
double cosTheta() const
int mod() const
double adc() const
double tdc() const
double rHit() const
int event() const
int qual() const
bool cutEtf()
Definition cut.cxx:308
void initial()
bool cutBarrel()
Definition cut.cxx:3
bool cutEndcap()
Definition cut.cxx:223
void setData(TTree *, unsigned int)
void setEtfDataFiles(std::vector< std::string > &)
void setEtfData(RecBTofCalHitCol &)
void setEndcapDataFiles(std::vector< std::string > &)
void setBarrelData(RecBTofCalHitCol &)
void setBarrelDataFiles(std::vector< std::string > &)
void setEndcapData(RecETofCalHitCol &)
double precision pisqo6 one
Definition qlconstants.h:4
Int_t run
Definition TofDataSet.h:18
Double_t tright
Definition TofDataSet.h:25
Int_t tofid
Definition TofDataSet.h:20
Double_t qright
Definition TofDataSet.h:23
Double_t path
Definition TofDataSet.h:29
Double_t texp
Definition TofDataSet.h:28
Double_t phi
Definition TofDataSet.h:30
Int_t event
Definition TofDataSet.h:19
Double_t dt
Definition TofDataSet.h:27
Int_t hitcase
Definition TofDataSet.h:35
Int_t strip
Definition TofDataSet.h:21
Double_t zrhit
Definition TofDataSet.h:26
Double_t tleft
Definition TofDataSet.h:24
Double_t theta
Definition TofDataSet.h:31
Double_t p
Definition TofDataSet.h:32
Double_t qleft
Definition TofDataSet.h:22