BOSS 7.0.6
BESIII Offline Software System
Loading...
Searching...
No Matches
calib_barrel_sigma.cxx
Go to the documentation of this file.
2#include "TF1.h"
3
4//the fit function to the single-end time resolution versus z
5/*
6static double singleEndFunc(double *x, double *par) {
7 double xx = x[0];
8 double f = par[0]*barrelRadius/sqrt(pow(barrelRadius,2)+pow(xx,2))+par[1]*xx+par[2]*pow(xx,2)+par[3]*pow(xx,3);
9 return f;
10}
11
12
13//the fit function to the double-end time resolution versus z
14static double doubleEndFunc(double *x, double *par)
15{
16 double xx = x[0];
17 double f = par[0]*pow(xx,2)/sqrt(pow(barrelRadius,2)+pow(xx,2))+par[1]+par[2]*pow(xx,2);
18 return f;
19}
20*/
21
23
24 nKind = 5; // 0:tleft, 1:tright, 2:t0, 3:tplus, 4:tminus
25 nBinPerCounter = nzbin;
26
29 CanvasPerCounterName.push_back( static_cast<string>("Barrel-offset") );
30 CanvasPerCounterName.push_back( static_cast<string>("Offset-TimeCorrelation") );
31 CanvasPerCounterName.push_back( static_cast<string>("Barrel-sigma") );
32 CanvasPerCounterName.push_back( static_cast<string>("Sigma-TimeCorrelation") );
33 nGraphPerCanvasPerCounter.push_back(3);
34 nGraphPerCanvasPerCounter.push_back(2);
35 nGraphPerCanvasPerCounter.push_back(3);
36 nGraphPerCanvasPerCounter.push_back(3);
37
38 nHistogram = 0;
39 nCanvas = 0;
40
41 int numGraphs = 0;
42 std::vector<unsigned int>::iterator iter = nGraphPerCanvasPerCounter.begin();
43 for( ; iter!=nGraphPerCanvasPerCounter.end(); iter++ ) {
44 numGraphs = numGraphs + (*iter);
45 }
46 if( numGraphs != nGraphTotalSigma ) {
47 cout << "tofcalgsec::calib_barrel_sigma: the number of Graphs is NOT reasonable!!!" << endl;
48 exit(0);
49 }
50
51 m_name = string("calib_barrel_sigma");
52
53 const int tbin = 150;
54 const double tbegin = -1.5;
55 const double tend = 1.5;
56
57 // histograms per counter
58 char hname[256];
59 for( unsigned int i=0; i<NBarrel; i++ ) {
60 m_result.push_back( HepVector(nBarrelSigma,0) );
61 for( unsigned int j=0; j<nKind; j++ ) {
62 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
63 if( j==0 ) { sprintf( hname, "tleft-id%i-z%i", i, k); }
64 else if( j==1 ) { sprintf( hname, "tright-id%i-z%i", i, k); }
65 else if( j==2 ) { sprintf( hname, "t0-id%i-z%i", i, k); }
66 else if( j==3 ) { sprintf( hname, "tplus-id%i-z%i", i, k); }
67 else if( j==4 ) { sprintf( hname, "tminus-id%i-z%i", i, k); }
68 m_histograms.push_back( new TH1F( hname, hname, tbin, tbegin, tend ) );
69
70 m_fitresult.push_back( HepVector(nParSigma,0) );
71 }
72 }
73 }
74
75 zpos.resize( nBinPerCounter );
76 zposerr.resize( nBinPerCounter );
77 zstep = ( zend - zbegin )/nBinPerCounter;
78 for( unsigned int i=0; i<nBinPerCounter; i++ ) {
79 zpos[i] = zbegin + ( i+0.5 )*zstep;
80 zposerr[i] = 0.5*zstep;
81 }
82
83}
84
85
87 m_fitresult.clear();
88 zpos.clear();
89 zposerr.clear();
90}
91
92
93void calib_barrel_sigma::calculate( RecordSet*& data, unsigned int icounter ) {
94
95 std::cout << setiosflags(ios::left) << setw(10) << icounter << setw(8) << data->size() << setw(30) << name() << std::endl;
96
97 if( data->size() > 0 ) {
98 std::vector<Record*>::iterator iter = data->begin();
99 for( ; iter!=data->end(); iter++ ) {
100 fillRecord( (*iter), icounter );
101 }
102 }
103 fitHistogram( icounter );
104 fillGraph( icounter );
105 fitGraph( icounter );
106
107 if( data->size() > 0 ) {
108 std::vector<Record*>::iterator iter = data->begin();
109 for( ; iter!=data->end(); iter++ ) {
110 updateData( (*iter), icounter );
111 fillRecordT0( (*iter), icounter );
112 }
113 }
114 fitHistogramT0( icounter );
115 fillGraphT0( icounter );
116 fitGraphT0( icounter );
117
118 return;
119}
120
121
122void calib_barrel_sigma::fillRecord( const Record* r, unsigned int icounter ) {
123
124 double zhit = r->zrhit();
125 if( zhit<zbegin || zhit>zend ) return;
126 int zbin = static_cast<int>((zhit-zbegin)/zstep);
127 if( zbin<0 ) { zbin = 0; }
128 else if( zbin>static_cast<int>(nBinPerCounter-1) ) {
129 cout << "tofcalgsec::calib_barrel_sigma:fillRecord: zhit is out of range, zhit=" << zhit << " zbin=" << zbin << endl;
130 return;
131 }
132
133 std::vector<TH1F*>::iterator iter = m_histograms.begin();
134 unsigned int number = icounter*nKind*nBinPerCounter + zbin;
135 (*(iter+number))->Fill( r->tleft() );
136 (*(iter+nBinPerCounter+number))->Fill( r->tright() );
137 (*(iter+3*nBinPerCounter+number))->Fill( (r->tleft()+r->tright())/2.0 );
138 (*(iter+4*nBinPerCounter+number))->Fill( (r->tleft()-r->tright())/2.0 );
139
140 return;
141}
142
143
144void calib_barrel_sigma::fitHistogram( unsigned int icounter ) {
145 TF1* g = new TF1("g", "gaus");
146 g->SetLineColor(2);
147 g->SetLineWidth(1);
148
149 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + icounter*nKind*nBinPerCounter;
150 std::vector<HepVector>::iterator iter2 = m_fitresult.begin() + icounter*nKind*nBinPerCounter;
151 for( unsigned int i=0; i<nKind; i++ ) {
152 for( unsigned int j=0; j<nBinPerCounter; j++ ) {
153 if( i!=2 ) {
154 (*iter1)->Fit( g, "Q");
155 (*iter2)[0] = g->GetParameter(1);
156 (*iter2)[1] = g->GetParError(1);
157 (*iter2)[2] = g->GetParameter(2);
158 (*iter2)[3] = g->GetParError(2);
159 }
160 iter1++;
161 iter2++;
162 }
163 }
164
165 return;
166
167}
168
169
170void calib_barrel_sigma::fillGraph( unsigned int icounter ) {
171
172 char gname1[256], gname2[256];
173
174 // fill graphs
175 // 4 canvas per counter,
176 // 1. offset of tleft, tright and t0 vs z
177 // 2. sigma of tleft, tright and t0 vs z
178 // 3. offset of tplus and tminus vs z
179 // 4. sigma of tplus, tminus and T_Correlation vs z
180 std::vector<double> toffset, toffseterr;
181 std::vector<double> tsigma, tsigmaerr;
182 toffset.resize( nBinPerCounter );
183 toffseterr.resize( nBinPerCounter );
184 tsigma.resize( nBinPerCounter );
185 tsigmaerr.resize( nBinPerCounter );
186
187 unsigned int number = 0;
188 std::vector<HepVector>::iterator iter = m_fitresult.begin() + icounter*nKind*nBinPerCounter;
189 for( unsigned int j=0; j<nKind; j++ ) {
190 if( j==0 ) { sprintf( gname1, "tleft-offset-tofid-%i", icounter ); }
191 else if( j==1 ) { sprintf( gname1, "tright-offset-tofid-%i", icounter ); }
192 else if( j==2 ) { sprintf( gname1, "t0-offset-tofid-%i", icounter ); }
193 else if( j==3 ) { sprintf( gname1, "tplus-offset-tofid-%i", icounter ); }
194 else if( j==4 ) { sprintf( gname1, "tminus-offset-tofid-%i", icounter ); }
195
196 m_graphs.push_back( new TH1F( gname1, gname1, nBinPerCounter, zbegin, zend ) );
197 std::vector<TH1F*>::iterator itgraph = m_graphs.end() - 1;
198
199 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
200 number = j*nBinPerCounter + k;
201 toffset[k] = (*(iter+number))[0];
202 toffseterr[k] = (*(iter+number))[1];
203 (*itgraph)->SetBinContent( k+1, toffset[k] );
204 (*itgraph)->SetBinError( k+1, toffseterr[k] );
205 }
206 (*itgraph)->SetMarkerSize(1.5);
207 if( j==0 || j==3) {
208 (*itgraph)->SetMarkerStyle(20);
209 (*itgraph)->SetMarkerColor(2);
210 (*itgraph)->SetMaximum( 0.15 );
211 (*itgraph)->SetMinimum(-0.15 );
212 }
213 else if( j==1 || j==4 ) {
214 (*itgraph)->SetMarkerStyle(21);
215 (*itgraph)->SetMarkerColor(4);
216 }
217 else if( j==2 ) {
218 (*itgraph)->SetMarkerStyle(4);
219 (*itgraph)->SetMarkerColor(2);
220 }
221 }
222
223 for( unsigned int j=0; j<nKind; j++ ) {
224 if( j==0 ) { sprintf( gname2, "tleft-sigma-tofid-%i", icounter ); }
225 else if( j==1 ) { sprintf( gname2, "tright-sigma-tofid-%i", icounter ); }
226 else if( j==2 ) { sprintf( gname2, "t0-sigma-tofid-%i", icounter ); }
227 else if( j==3 ) { sprintf( gname2, "tplus-sigma-tofid-%i", icounter ); }
228 else if( j==4 ) { sprintf( gname2, "tminus-sigma-tofid-%i", icounter ); }
229 m_graphs.push_back( new TH1F( gname2, gname2, nBinPerCounter, zbegin, zend ) );
230 std::vector<TH1F*>::iterator itgraph = m_graphs.end() - 1;
231
232 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
233 number = j*nBinPerCounter + k;
234 tsigma[k] = (*(iter+number))[2];
235 tsigmaerr[k] = (*(iter+number))[3];
236 (*itgraph)->SetBinContent( k+1, tsigma[k] );
237 (*itgraph)->SetBinError( k+1, tsigmaerr[k] );
238 }
239 (*itgraph)->SetMarkerSize(1.5);
240 if( j==0 || j==3 ) {
241 (*itgraph)->SetMarkerStyle(20);
242 (*itgraph)->SetMarkerColor(2);
243 (*itgraph)->SetMaximum( 0.3 );
244 (*itgraph)->SetMinimum( 0.0 );
245 }
246 else if( j==1 || j==4 ) {
247 (*itgraph)->SetMarkerStyle(21);
248 (*itgraph)->SetMarkerColor(4);
249 }
250 else if( j==2 ) {
251 (*itgraph)->SetMarkerStyle(4);
252 (*itgraph)->SetMarkerColor(2);
253 }
254 }
255
256 sprintf( gname2, "sigma-tofid-%i", icounter );
257 m_graphs.push_back( new TH1F( gname2, gname2, nBinPerCounter, zbegin, zend ) );
258 std::vector<TH1F*>::iterator itgraph = m_graphs.end() - 1;
259 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
260 number = (nKind-1)*nBinPerCounter + k;
261 double sigPlus = (*(iter+number-nBinPerCounter))[2];
262 double sigMinus = (*(iter+number))[2];
263 double sigErrPlus = (*(iter+number-nBinPerCounter))[3];
264 double sigErrMinus = (*(iter+number))[3];
265 tsigma[k] = sqrt( sigPlus*sigPlus - sigMinus*sigMinus );
266 tsigmaerr[k] = sqrt( sigErrPlus*sigErrPlus + sigErrMinus*sigErrMinus );
267 (*itgraph)->SetBinContent( k+1, tsigma[k] );
268 (*itgraph)->SetBinError( k+1, tsigmaerr[k] );
269 }
270 (*itgraph)->SetMarkerSize(1.5);
271 (*itgraph)->SetMarkerStyle(4);
272 (*itgraph)->SetMarkerColor(2);
273
274 return;
275}
276
277
278void calib_barrel_sigma::fitGraph( unsigned int icounter ) {
279
280 TF1* fsingle = new TF1("fsingle", "pol4");
281 fsingle->SetLineColor(1);
282 fsingle->SetLineWidth(1);
283
284 std::vector<unsigned int>::iterator itnumber = nGraphPerCanvasPerCounter.begin();
285 std::vector<TH1F*>::iterator itgraph = m_graphs.begin() + icounter*nGraphTotalSigma + (*itnumber) + (*(itnumber+1));
286
287 fsingle->SetParameter( 0, 0.14 );
288 fsingle->SetParameter( 1, -4.0e-4 );
289 (*itgraph)->Fit( "fsingle", "QR", "", zbegin, zend );
290 X = HepVector( m_npar, 0 );
291 for( unsigned int i=0; i<5; i++ ) {
292 X[i] = fsingle->GetParameter(i);
293 }
294
295 fsingle->SetParameter( 0, 0.14 );
296 fsingle->SetParameter( 1, 4.0e-4 );
297 (*(itgraph+1))->Fit( "fsingle", "QR", "", zbegin, zend );
298 for( unsigned int i=0; i<5; i++ ) {
299 X[i+5] = fsingle->GetParameter(i);
300 }
301
302 std::vector<HepVector>::iterator iter = m_result.begin() + icounter;
303 (*iter) = X;
304
305 return;
306}
307
308
309void calib_barrel_sigma::updateData( Record* r, unsigned int icounter ) {
310 double zhit = r->zrhit();
311 double t1 = r->tleft();
312 double t2 = r->tright();
313
314 double par1[5], par2[5];
315 std::vector<HepVector>::iterator iter = m_result.begin() + icounter;
316 for( unsigned int i=0; i<5; i++ ) {
317 par1[i] = (*iter)[i];
318 par2[i] = (*iter)[i+5];
319 }
320
321 double tsigma1 = par1[0]+par1[1]*zhit+par1[2]*pow(zhit,2)+par1[3]*pow(zhit,3) + par1[4]*pow(zhit,4);
322 double tsigma2 = par2[0]+par2[1]*zhit+par2[2]*pow(zhit,2)+par2[3]*pow(zhit,3) + par2[4]*pow(zhit,4);
323 double tc = m_tcorrelation[0];
324
325 double weight1 = (tsigma2*tsigma2-tc*tc)/(tsigma1*tsigma1+tsigma2*tsigma2-2.0*tc*tc);
326 double weight2 = (tsigma1*tsigma1-tc*tc)/(tsigma1*tsigma1+tsigma2*tsigma2-2.0*tc*tc);
327 double t0 = weight1*t1 + weight2*t2;
328
329 r->setT0( t0 );
330
331 return;
332}
333
334
335void calib_barrel_sigma::fillRecordT0( const Record* r, unsigned int icounter ) {
336 double zhit = r->zrhit();
337 if( zhit<zbegin || zhit>zend ) return;
338 int zbin = static_cast<int>((zhit-zbegin)/zstep);
339 if( zbin<0 ) { zbin = 0; }
340 else if( zbin>static_cast<int>(nBinPerCounter-1) ) {
341 cout << "tofcalgsec::calib_barrel_sigma:fillRecordT0: zhit is out of range, zhit=" << zhit << " zbin=" << zbin << endl;
342 return;
343 }
344
345 std::vector<TH1F*>::iterator iter = m_histograms.begin();
346 unsigned int number = icounter*nKind*nBinPerCounter + 2*nBinPerCounter + zbin;
347 (*(iter+number))->Fill( r->t0() );
348
349 return;
350}
351
352
353void calib_barrel_sigma::fitHistogramT0( unsigned int icounter ) {
354 TF1* g = new TF1("g", "gaus");
355 g->SetLineColor(2);
356 g->SetLineWidth(1);
357
358 std::vector<TH1F*>::iterator iter1 = m_histograms.begin() + icounter*nKind*nBinPerCounter + 2*nBinPerCounter;
359 std::vector<HepVector>::iterator iter2 = m_fitresult.begin() + icounter*nKind*nBinPerCounter + 2*nBinPerCounter;
360 for( unsigned int j=0; j<nBinPerCounter; j++, iter1++, iter2++ ) {
361 (*iter1)->Fit( g, "Q");
362 (*iter2)[0] = g->GetParameter(1);
363 (*iter2)[1] = g->GetParError(1);
364 (*iter2)[2] = g->GetParameter(2);
365 (*iter2)[3] = g->GetParError(2);
366 }
367
368 return;
369}
370
371
372void calib_barrel_sigma::fillGraphT0( unsigned int icounter ) {
373 char gname1[256], gname2[256];
374
375 // fill graphs
376 // 4 canvas per counter,
377 // 1. offset of tleft, tright and t0 vs z
378 // 2. sigma of tleft, tright and t0 vs z
379 // 3. offset of tplus and tminus vs z
380 // 4. sigma of tplus, tminus and T_Correlation vs z
381 std::vector<double> toffset, toffseterr;
382 std::vector<double> tsigma, tsigmaerr;
383 toffset.resize( nBinPerCounter );
384 toffseterr.resize( nBinPerCounter );
385 tsigma.resize( nBinPerCounter );
386 tsigmaerr.resize( nBinPerCounter );
387
388 std::vector<HepVector>::iterator iter = m_fitresult.begin() + icounter*nKind*nBinPerCounter + 2*nBinPerCounter;
389 for( unsigned int k=0; k<nBinPerCounter; k++ ) {
390 toffset[k] = (*(iter+k))[0];
391 toffseterr[k] = (*(iter+k))[1];
392 tsigma[k] = (*(iter+k))[2];
393 tsigmaerr[k] = (*(iter+k))[3];
394 }
395
396 sprintf( gname1, "offset-tofid-%i", icounter );
397 std::vector<TH1F*>::iterator itgraph = m_graphs.begin() + icounter*nGraphTotalSigma + 2;
398 for( unsigned int i=0; i<nBinPerCounter; i++ ) {
399 // (*itgraph)->SetPoint( i, zpos[i], toffset[i] );
400 // (*itgraph)->SetPointError( i, zposerr[i], toffseterr[i] );
401 (*itgraph)->SetBinContent( i+1, toffset[i] );
402 (*itgraph)->SetBinError( i+1, toffseterr[i] );
403 }
404
405 sprintf( gname2, "sigma-tofid-%i", icounter );
406 itgraph = m_graphs.begin() + icounter*nGraphTotalSigma + 7;
407 for( unsigned int i=0; i<nBinPerCounter; i++ ) {
408 // (*itgraph)->SetPoint( i, zpos[i], tsigma[i] );
409 // (*itgraph)->SetPointError( i, zposerr[i], tsigmaerr[i] );
410 (*itgraph)->SetBinContent( i+1, tsigma[i] );
411 (*itgraph)->SetBinError( i+1, tsigmaerr[i] );
412 }
413
414 return;
415}
416
417
418void calib_barrel_sigma::fitGraphT0( unsigned int icounter ) {
419
420 // TF1 *fdouble = new TF1( "fdouble", doubleEndFunc, zbegin, zend, 3 );
421 TF1 *fdouble = new TF1( "fdouble", "pol4", zbegin, zend );
422 fdouble->SetLineColor(1);
423 fdouble->SetLineWidth(1);
424
425 std::vector<unsigned int>::iterator itnumber = nGraphPerCanvasPerCounter.begin();
426 std::vector<TH1F*>::iterator itgraph = m_graphs.begin() + icounter*nGraphTotalSigma + (*itnumber) + (*(itnumber+1)) + 2;
427 (*itgraph)->Fit( "fdouble", "Q", "", zbegin, zend );
428
429 std::vector<HepVector>::iterator iter = m_result.begin() + icounter;
430 (*iter)[10] = fdouble->GetParameter(0);
431 (*iter)[11] = fdouble->GetParameter(1);
432 (*iter)[12] = fdouble->GetParameter(2);
433 (*iter)[13] = fdouble->GetParameter(3);
434 (*iter)[14] = fdouble->GetParameter(4);
435
436 return;
437}
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
const double zend
Definition: TofCalibFit.h:13
const double zbegin
Definition: TofCalibFit.h:12
std::vector< Record * > RecordSet
Definition: TofDataSet.h:98
const unsigned int NBarrel
Definition: TofDataSet.h:12
const int nParSigma
const int nBarrelSigma
const int nGraphTotalSigma
double tleft() const
Definition: TofDataSet.h:59
double tright() const
Definition: TofDataSet.h:60
void setT0(double t0)
Definition: TofDataSet.h:75
double t0() const
Definition: TofDataSet.h:68
double zrhit() const
Definition: TofDataSet.h:61
string m_name
Definition: TofCalibFit.h:50
std::vector< HepVector > m_result
Definition: TofCalibFit.h:55
unsigned int nCanvas
Definition: TofCalibFit.h:47
std::vector< TH1F * > m_histograms
Definition: TofCalibFit.h:53
unsigned int nBinPerCounter
Definition: TofCalibFit.h:41
HepVector m_tcorrelation
Definition: TofCalibFit.h:60
unsigned int nKind
Definition: TofCalibFit.h:40
std::vector< unsigned int > nGraphPerCanvasPerCounter
Definition: TofCalibFit.h:45
unsigned int nCanvasPerCounter
Definition: TofCalibFit.h:44
unsigned int nHistPerCounter
Definition: TofCalibFit.h:43
const string & name() const
Definition: TofCalibFit.h:27
HepVector X
Definition: TofCalibFit.h:51
unsigned int nHistogram
Definition: TofCalibFit.h:46
std::vector< TH1F * > m_graphs
Definition: TofCalibFit.h:54
std::vector< string > CanvasPerCounterName
Definition: TofCalibFit.h:57
calib_barrel_sigma(const unsigned int nzbin)
void calculate(RecordSet *&data, unsigned int icounter)
sprintf(cut,"kal_costheta0_em>-0.93&&kal_costheta0_em<0.93&&kal_pxy0_em>=0.05+%d*0.1&&kal_pxy0_em<0.15+%d*0.1&&NGch>=2", j, j)