CGEM BOSS 6.6.5.h
BESIII Offline Software System
Loading...
Searching...
No Matches
TofDataSet.cxx
Go to the documentation of this file.
2
4 initial();
5 return;
6}
7
8
10
11
13 m_run = 0;
14 m_event = -1;
15 m_tofid = -1;
16 m_qleft = -99.0;
17 m_qright = -99.0;
18 m_tleft = -99.0;
19 m_tright = -99.0;
20 m_zrhit = -999.0;
21 m_texp = -99.0;
22 m_path = -99.0;
23 m_phi = -99.0;
24 m_theta = -99.0;
25 m_p = -9.0;
26 m_t0 = -99.0;
27 m_q0 = -99.0;
28 m_hitcase = -1;
29 return;
30}
31
32
34 initial();
35 if( hit ) {
36 m_run = hit->run();
37 m_event = hit->event();
38 m_tofid = hit->mod();
39 m_qleft = hit->adc1();
40 m_qright = hit->adc2();
41 m_tleft = hit->tdc1();
42 m_tright = hit->tdc2();
43 m_zrhit = hit->zHit();
44 m_texp = hit->tpred();
45 m_path = hit->path();
46 m_phi = hit->deltaPhi();
47 m_theta = hit->sinTheta();
48 m_p = hit->p();
49 m_hitcase = hit->qual();
50 }
51 return;
52}
53
54
56 initial();
57 if( hit ) {
58 m_run = hit->run();
59 m_event = hit->event();
60 m_tofid = hit->mod();
61 m_qleft = hit->adc();
62 m_tleft = hit->tdc();
63 m_zrhit = hit->rHit();
64 m_texp = hit->tpred();
65 m_path = hit->path();
66 m_phi = hit->deltaPhi();
67 m_theta = hit->cosTheta();
68 m_p = hit->p();
69 m_hitcase = hit->qual();
70 }
71 return;
72}
73
74
76 initial();
77 if( &one ) {
78 m_run = one.run;
79 m_event = one.event;
80 m_tofid = one.tofid;
81 m_qleft = one.qleft;
82 m_tleft = one.tleft;
83 if( one.hitcase >= 0 && one.hitcase <=2 ) {
84 m_qright = one.qright;
85 m_tright = one.tright;
86 }
87 m_zrhit = one.zrhit;
88 m_texp = one.texp;
89 m_path = one.path;
90 m_phi = one.phi;
91 m_theta = one.theta;
92 m_p = one.p;
93 m_hitcase = one.hitcase;
94 }
95 return;
96}
97
98
99void TofDataSet::setData( TTree* t, bool isBarrel ) {
100 if( t && t->GetEntries()>0 ) {
101 rootRecord item;
102 t->SetBranchAddress( "run", &item.run );
103 t->SetBranchAddress( "event", &item.event );
104 t->SetBranchAddress( "tofid", &item.tofid );
105 t->SetBranchAddress( "qleft", &item.qleft );
106 t->SetBranchAddress( "qright", &item.qright );
107 t->SetBranchAddress( "tleft", &item.tleft );
108 t->SetBranchAddress( "tright", &item.tright );
109 t->SetBranchAddress( "zrhit", &item.zrhit );
110 t->SetBranchAddress( "texp", &item.texp );
111 t->SetBranchAddress( "path", &item.path );
112 t->SetBranchAddress( "phi", &item.phi );
113 t->SetBranchAddress( "theta", &item.theta );
114 t->SetBranchAddress( "p", &item.p );
115 t->SetBranchAddress( "hitcase", &item.hitcase );
116
117 for( unsigned int i=0; i<t->GetEntries(); i++ ) {
118 t->GetEntry(i);
119 if( isBarrel && item.hitcase>=0 && item.hitcase<=2 ) {
120 Record *r = new Record( item );
121 if( r->cutBarrel() ) {
122 unsigned int tofID = item.tofid;
123 barrelData[tofID]->push_back(r);
124 }
125 else {
126 delete r;
127 }
128 }
129 else if( !isBarrel && ( item.hitcase==3 || item.hitcase==4 ) ) {
130 Record *r = new Record( item );
131 if( r->cutEndcap() ) {
132 unsigned int tofID = item.tofid;
133 endcapData[tofID]->push_back(r);
134 }
135 else {
136 delete r;
137 }
138 }
139 }
140 }
141 else {
142 std::cerr << "Error: a invalid tree or a blank tree, When converting a tree to TofDataSet,exit" << std::endl;
143 }
144
145 return;
146}
147
148
150 for( unsigned int i=0; i<NBarrel; i++ ) {
151 barrelData[i] = new RecordSet;
152 }
153 for( unsigned int i=0; i<NEndcap; i++ ) {
154 endcapData[i] = new RecordSet;
155 }
156 return;
157}
158
159
161 for( unsigned int i=0; i<NBarrel; i++ ) {
162 barrelData[i]->clear();
163 delete barrelData[i];
164 }
165 for( unsigned int i=0; i<NEndcap; i++ ) {
166 endcapData[i]->clear();
167 delete endcapData[i];
168 }
169 return;
170}
171
172
173void TofDataSet::setBarrelDataFiles( std::vector<std::string>& barrelFiles ) {
174 TChain* data_barrel = new TChain("btrk");
175 if( !data_barrel ) {
176 std::cerr << " tofcalgsec Error Msg: creating a tree[barrel] fails in TofDataSet()"<<std::endl;
177 throw "Error Msg: creating a tree fails in TofDataSet() ";
178 }
179 std::cout<<"begin reading barrel data file ... "<<std::endl;
180 try{
181 for( std::vector<std::string>::iterator it=barrelFiles.begin(); it!=barrelFiles.end(); it++ ) {
182 std::cout << " Add file : " << (*it) << std::endl;
183 data_barrel->Add( (*it).c_str() );
184 }
185 }
186 catch(...){
187 std::cerr << "tofcalgsec Error Msg : in TofDataSet::setDataFiles(std::vector<std::string>&) " << std::endl;
188 return;
189 }
190 setData( data_barrel, true );
191 delete data_barrel;
192 return;
193}
194
195
196void TofDataSet::setEndcapDataFiles( std::vector<std::string>& endcapFiles ) {
197 TChain* data_endcap = new TChain("etrk");
198 if( !data_endcap ) {
199 std::cerr << " tofcalgsec Error Msg: creating a tree[endcap] fails in TofDataSet()"<<std::endl;
200 throw "Error Msg: creating a tree fails in TofDataSet() ";
201 }
202 std::cout<<"begin reading endcap data file ... "<<std::endl;
203 try{
204 for( std::vector<std::string>::iterator it=endcapFiles.begin(); it!=endcapFiles.end(); it++ ) {
205 std::cout << " Add file : " << (*it) << std::endl;
206 data_endcap->Add( (*it).c_str() );
207 }
208 }
209 catch(...){
210 std::cerr << "tofcalgsec Error Msg : in TofDataSet::setDataFiles(std::vector<std::string>&) " << std::endl;
211 return;
212 }
213 setData( data_endcap, false );
214 delete data_endcap;
215 return;
216}
217
218
220 RecBTofCalHitCol::iterator iter = bhitcol.begin();
221 for( ; iter!=bhitcol.end(); iter++ ) {
222 int tofid = (*iter)->mod();
223 if( tofid<0 || tofid>175 ) continue;
224
225 Record *r = new Record( (*iter) );
226 if( r->cutBarrel() ) {
227 barrelData[tofid]->push_back(r);
228 }
229 else {
230 delete r;
231 }
232 }
233 return;
234}
235
236
238 RecETofCalHitCol::iterator iter = ehitcol.begin();
239 for( ; iter!=ehitcol.end(); iter++ ) {
240 int tofid = (*iter)->mod();
241 if( tofid<0 || tofid>95 ) continue;
242
243 Record *r = new Record( (*iter) );
244 if( r->cutEndcap() ) {
245 endcapData[tofid]->push_back(r);
246 }
247 else {
248 delete r;
249 }
250 }
251 return;
252}
const DifNumber one
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
ObjectVector< RecBTofCalHit > RecBTofCalHitCol
ObjectVector< RecETofCalHit > RecETofCalHitCol
std::vector< Record * > RecordSet
Definition TofDataSet.h:89
const unsigned int NBarrel
Definition TofDataSet.h:12
const unsigned int NEndcap
Definition TofDataSet.h:13
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
void initial()
bool cutBarrel()
Definition cut.cxx:3
bool cutEndcap()
Definition cut.cxx:88
void setEndcapDataFiles(std::vector< std::string > &)
void setBarrelData(RecBTofCalHitCol &)
void setData(TTree *, bool)
void setBarrelDataFiles(std::vector< std::string > &)
void setEndcapData(RecETofCalHitCol &)
Int_t run
Definition TofDataSet.h:16
Double_t tright
Definition TofDataSet.h:22
Int_t tofid
Definition TofDataSet.h:18
Double_t qright
Definition TofDataSet.h:20
Double_t path
Definition TofDataSet.h:25
Double_t texp
Definition TofDataSet.h:24
Double_t phi
Definition TofDataSet.h:26
Int_t event
Definition TofDataSet.h:17
Int_t hitcase
Definition TofDataSet.h:31
Double_t zrhit
Definition TofDataSet.h:23
Double_t tleft
Definition TofDataSet.h:21
Double_t theta
Definition TofDataSet.h:27
Double_t p
Definition TofDataSet.h:28
Double_t qleft
Definition TofDataSet.h:19
int t()
Definition t.c:1