BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
XtInteCalib.cpp
Go to the documentation of this file.
2#include "include/fun.h"
3
4#include <cmath>
5
6#include "TF1.h"
7#include "TMinuit.h"
8#include "TTree.h"
9
11 m_nMaxGrPoint = 5000;
12 for ( int lay = 0; lay < NLAYER; lay++ )
13 {
14 m_tbinWid[lay][0] = 5.0;
15 m_tbinWid[lay][1] = 10.0;
16 m_tbinWid[lay][2] = 20.0;
17
18 m_tbinLim[lay][0] = -10.0;
19 m_tbinLim[lay][1] = 30.0;
20 if ( lay < 8 ) m_tbinLim[lay][2] = 210.0;
21 else m_tbinLim[lay][2] = 350.0;
22 m_tbinLim[lay][3] = 990.0;
23 }
24 cout << "Calibration type: XtInteCalib" << endl;
25}
26
28
29void XtInteCalib::init( TObjArray* hlist, MdcCosGeom* pGeom ) {
30 CalibBase::init( hlist, pGeom );
31
32 m_fdPf = new TFolder( "mfdProfile", "fdProfile" );
33 hlist->Add( m_fdPf );
34
35 m_haxis = new TH2F( "axis", "", 200, -50, 1000, 50, 0, 10 );
36 m_haxis->SetStats( 0 );
37 m_fdPf->Add( m_haxis );
38
39 char hname[200];
40 for ( int lay = 0; lay < NLAYER; lay++ )
41 {
42 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
43 {
44 for ( int lr = 0; lr < 2; lr++ )
45 {
46 sprintf( hname, "mxt%02d_%02d_%d_gr", lay, iEntr, lr );
47 m_grXt[lay][iEntr][lr] = new TGraph();
48 m_grXt[lay][iEntr][lr]->SetName( hname );
49 m_grXt[lay][iEntr][lr]->SetMarkerColor( 2 );
50 m_fdPf->Add( m_grXt[lay][iEntr][lr] );
51
52 int xbinN =
53 (int)( ( m_tbinLim[lay][1] - m_tbinLim[lay][0] ) / m_tbinWid[lay][0] + 0.5 );
54 sprintf( hname, "mxt%02d_%02d_%d_near", lay, iEntr, lr );
55 m_pfNear[lay][iEntr][lr] =
56 new TProfile( hname, hname, xbinN, m_tbinLim[lay][0], m_tbinLim[lay][1] );
57 m_fdPf->Add( m_pfNear[lay][iEntr][lr] );
58
59 int xbinM =
60 (int)( ( m_tbinLim[lay][2] - m_tbinLim[lay][1] ) / m_tbinWid[lay][1] + 0.5 );
61 sprintf( hname, "mxt%02d_%02d_%d_mid", lay, iEntr, lr );
62 m_pfMid[lay][iEntr][lr] =
63 new TProfile( hname, hname, xbinM, m_tbinLim[lay][1], m_tbinLim[lay][2] );
64 m_fdPf->Add( m_pfMid[lay][iEntr][lr] );
65
66 int xbinF =
67 (int)( ( m_tbinLim[lay][3] - m_tbinLim[lay][2] ) / m_tbinWid[lay][2] + 0.5 );
68 sprintf( hname, "mxt%02d_%02d_%d_far", lay, iEntr, lr );
69 m_pfFar[lay][iEntr][lr] =
70 new TProfile( hname, hname, xbinF, m_tbinLim[lay][2], m_tbinLim[lay][3] );
71 m_fdPf->Add( m_pfFar[lay][iEntr][lr] );
72 }
73 }
74 }
75}
76
77void XtInteCalib::mergeHist( TFile* fhist ) {
78 CalibBase::mergeHist( fhist );
79
80 double tdr, doca;
81 char hname[200];
82 TProfile* pr;
83 TFolder* fd = (TFolder*)fhist->Get( "fdProfile" );
84 for ( int lay = 0; lay < NLAYER; lay++ )
85 {
86 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
87 {
88 for ( int lr = 0; lr < 2; lr++ )
89 {
90 if ( ( m_grXt[lay][iEntr][lr]->GetN() ) < m_nMaxGrPoint )
91 {
92 sprintf( hname, "xt%02d_%02d_%d_gr", lay, iEntr, lr );
93 TGraph* gr = (TGraph*)fd->FindObjectAny( hname );
94 int nPoint = gr->GetN();
95 for ( int i = 0; i < nPoint; i++ )
96 {
97 gr->GetPoint( i, tdr, doca );
98 int np = m_grXt[lay][iEntr][lr]->GetN();
99 m_grXt[lay][iEntr][lr]->SetPoint( np, tdr, doca );
100 }
101 }
102 sprintf( hname, "xt%02d_%02d_%d_near", lay, iEntr, lr );
103 pr = (TProfile*)fd->FindObjectAny( hname );
104 m_pfNear[lay][iEntr][lr]->Add( pr );
105
106 sprintf( hname, "xt%02d_%02d_%d_mid", lay, iEntr, lr );
107 pr = (TProfile*)fd->FindObjectAny( hname );
108 m_pfMid[lay][iEntr][lr]->Add( pr );
109
110 sprintf( hname, "xt%02d_%02d_%d_far", lay, iEntr, lr );
111 pr = (TProfile*)fd->FindObjectAny( hname );
112 m_pfFar[lay][iEntr][lr]->Add( pr );
113 }
114 }
115 }
116}
117
118void XtInteCalib::calib( MdcCalibConst* calconst, TObjArray* newXtList, TObjArray* r2tList ) {
119 CalibBase::calib( calconst, newXtList, r2tList );
120 bool fgOldXt = saveOldXt( newXtList );
121 newXtList->Clear();
122
123 TGraph* grXtfit[NLAYER][NENTRXT][2];
124 char hname[200];
125 for ( int lay = 0; lay < NLAYER; lay++ )
126 {
127 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
128 {
129 for ( int lr = 0; lr < 2; lr++ )
130 {
131 m_vt.clear();
132 m_vd.clear();
133 m_entries.clear();
134 for ( int iPf = 0; iPf < 3; iPf++ )
135 {
136 TProfile* pro;
137 if ( 0 == iPf ) pro = m_pfNear[lay][iEntr][lr];
138 else if ( 1 == iPf ) pro = m_pfMid[lay][iEntr][lr];
139 else pro = m_pfFar[lay][iEntr][lr];
140
141 int nbin = pro->GetNbinsX();
142 for ( int i = 0; i < nbin; i++ )
143 {
144 double tt = pro->GetBinCenter( i + 1 );
145 double dd = pro->GetBinContent( i + 1 );
146 double entries = pro->GetBinEntries( i + 1 );
147 if ( entries > 10 )
148 {
149 m_vt.push_back( tt );
150 m_vd.push_back( dd );
151 m_entries.push_back( entries );
152 }
153 }
154 }
155 unsigned vsize = m_vt.size();
156 if ( vsize > 10 )
157 {
158 for ( int i = 0; i < 2; i++ )
159 {
160 double slope =
161 ( m_vd[vsize - 1] - m_vd[vsize - 2] ) / ( m_vt[vsize - 1] - m_vt[vsize - 2] );
162 if ( fabs( slope ) > 0.04 )
163 { // 0.8mm/20ns
164 m_vt.pop_back();
165 m_vd.pop_back();
166 m_entries.pop_back();
167 vsize = m_vt.size();
168 }
169 }
170 }
171 sprintf( hname, "grXtFit%02d_%02d_%d", lay, iEntr, lr );
172 grXtfit[lay][iEntr][lr] = new TGraph();
173 grXtfit[lay][iEntr][lr]->SetName( hname );
174 grXtfit[lay][iEntr][lr]->SetMarkerStyle( 20 );
175 m_fgFit[lay][iEntr][lr] = getXt( lay, iEntr, lr, grXtfit[lay][iEntr][lr] );
176 }
177 }
178 }
179
180 double tdr, doca;
181 for ( int lay = 0; lay < NLAYER; lay++ )
182 {
183 double tCut = 500.0;
184 if ( lay < 8 ) tCut = 400.0;
185 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
186 {
187 for ( int lr = 0; lr < 2; lr++ )
188 {
189 if ( !m_fgFit[lay][iEntr][lr] )
190 {
191 int iEntrNew = findXtEntr( lay, iEntr, lr );
192 if ( -1 != iEntrNew )
193 {
194 int npoint = grXtfit[lay][iEntrNew][lr]->GetN();
195 for ( int i = 0; i < npoint; i++ )
196 {
197 grXtfit[lay][iEntrNew][lr]->GetPoint( i, tdr, doca );
198 grXtfit[lay][iEntr][lr]->SetPoint( i, tdr, doca );
199 }
200 }
201 else if ( fgOldXt )
202 {
203 cout << grXtfit[lay][iEntr][lr]->GetName() << " use old x-t" << endl;
204 int npoint = m_grXtOld[lay][iEntr][lr]->GetN();
205 for ( int i = 0; i < npoint; i++ )
206 {
207 m_grXtOld[lay][iEntr][lr]->GetPoint( i, tdr, doca );
208 grXtfit[lay][iEntr][lr]->SetPoint( i, tdr, doca );
209 }
210 }
211 }
212 int nn = grXtfit[lay][iEntr][lr]->GetN();
213 double tmax, dmax;
214 grXtfit[lay][iEntr][lr]->GetPoint( nn - 1, tmax, dmax );
215 if ( tmax > tCut ) m_fgEdge[lay][iEntr][lr] = true;
216 else m_fgEdge[lay][iEntr][lr] = false;
217 }
218 }
219 }
220
221 for ( int lay = 0; lay < NLAYER; lay++ )
222 {
223 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
224 {
225 for ( int lr = 0; lr < 2; lr++ )
226 {
227 if ( !m_fgEdge[lay][iEntr][lr] )
228 {
229 int iEntrNew = findXtEntrEdge( lay, iEntr, lr );
230 if ( -1 != iEntrNew )
231 {
232 double t1, d1;
233 int n1 = grXtfit[lay][iEntr][lr]->GetN();
234 grXtfit[lay][iEntr][lr]->GetPoint( n1 - 1, t1, d1 );
235 double t2, d2;
236 int n2 = grXtfit[lay][iEntrNew][lr]->GetN();
237 for ( int i = 0; i < n2; i++ )
238 {
239 grXtfit[lay][iEntrNew][lr]->GetPoint( i, t2, d2 );
240 if ( t2 > t1 )
241 {
242 grXtfit[lay][iEntr][lr]->SetPoint( n1, t2, d2 );
243 n1++;
244 }
245 }
246 }
247 }
248 }
249 }
250 }
251
252 TTree* xttr[NLAYER][NENTRXT][2];
253 for ( int lay = 0; lay < NLAYER; lay++ )
254 {
255 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
256 {
257 for ( int lr = 0; lr < 2; lr++ )
258 {
259 sprintf( hname, "trNewXt%02d_%02d_%d", lay, iEntr, lr );
260 xttr[lay][iEntr][lr] = new TTree( hname, hname );
261 xttr[lay][iEntr][lr]->Branch( "t", &tdr, "t/D" );
262 xttr[lay][iEntr][lr]->Branch( "d", &doca, "d/D" );
263 if ( 0 == gFgCalib[lay] )
264 {
265 int npoint = m_grXtOld[lay][iEntr][lr]->GetN();
266 for ( int i = 0; i < npoint; i++ )
267 {
268 m_grXtOld[lay][iEntr][lr]->GetPoint( i, tdr, doca );
269 xttr[lay][iEntr][lr]->Fill();
270 }
271 }
272 else
273 {
274 int npoint = grXtfit[lay][iEntr][lr]->GetN();
275 for ( int i = 0; i < npoint; i++ )
276 {
277 grXtfit[lay][iEntr][lr]->GetPoint( i, tdr, doca );
278 xttr[lay][iEntr][lr]->Fill();
279 }
280 }
281 newXtList->Add( xttr[lay][iEntr][lr] );
282 }
283 }
284 }
285
286 for ( int lay = 0; lay < NLAYER; lay++ )
287 {
288 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
289 {
290 for ( int lr = 0; lr < 2; lr++ ) { delete grXtfit[lay][iEntr][lr]; }
291 }
292 }
293
294 // fnewXt.Close();
295 renameHist();
296}
297
298void XtInteCalib::renameHist() {
299 char hname[200];
300 m_fdPf->SetName( "fdProfile" );
301 for ( int lay = 0; lay < NLAYER; lay++ )
302 {
303 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
304 {
305 for ( int lr = 0; lr < 2; lr++ )
306 {
307 sprintf( hname, "xt%02d_%02d_%d_gr", lay, iEntr, lr );
308 m_grXt[lay][iEntr][lr]->SetName( hname );
309 sprintf( hname, "xt%02d_%02d_%d_near", lay, iEntr, lr );
310 m_pfNear[lay][iEntr][lr]->SetName( hname );
311 sprintf( hname, "xt%02d_%02d_%d_mid", lay, iEntr, lr );
312 m_pfMid[lay][iEntr][lr]->SetName( hname );
313 sprintf( hname, "xt%02d_%02d_%d_far", lay, iEntr, lr );
314 m_pfFar[lay][iEntr][lr]->SetName( hname );
315 }
316 }
317 }
318}
319
320bool XtInteCalib::saveOldXt( TObjArray* newXtList ) {
321 char hname[200];
322 Int_t entries = newXtList->GetEntries();
323 cout << "entries of newXtList " << entries << endl;
324 if ( entries < 1548 )
325 {
326 cout << "can not get old x-t" << endl;
327 return false;
328 }
329 for ( int lay = 0; lay < NLAYER; lay++ )
330 {
331 for ( int iEntr = 0; iEntr < NENTRXT; iEntr++ )
332 {
333 for ( int lr = 0; lr < 2; lr++ )
334 {
335 sprintf( hname, "grXtOld%02d_%02d_%d", lay, iEntr, lr );
336 m_grXtOld[lay][iEntr][lr] = new TGraph();
337 m_grXtOld[lay][iEntr][lr]->SetName( hname );
338
339 double t, d;
340 sprintf( hname, "trNewXt%02d_%02d_%d", lay, iEntr, lr );
341 TTree* tr = (TTree*)( newXtList->FindObject( hname ) );
342 // cout << setw(15) << tr->GetName() << setw(15) << tr->GetEntries()
343 // << endl;
344 tr->SetBranchAddress( "t", &t );
345 tr->SetBranchAddress( "d", &d );
346 int nPoint = tr->GetEntries();
347 if ( nPoint < 20 )
348 {
349 cout << "can not get old x-t: " << hname << endl;
350 return false;
351 }
352 for ( int i = 0; i < nPoint; i++ )
353 {
354 tr->GetEntry( i );
355 m_grXtOld[lay][iEntr][lr]->SetPoint( i, t, d );
356 }
357 }
358 }
359 }
360 // for(int lay=0; lay<NLAYER; lay++){
361 // for(int iEntr=0; iEntr<NENTRXT; iEntr++){
362 // for(int lr=0; lr<2; lr++){
363 // sprintf(hname, "trNewXt%02d_%02d_%d", lay, iEntr, lr);
364 // delete (newXtList->FindObject(hname));
365 // }
366 // }
367 // }
368 return true;
369}
370
371bool XtInteCalib::getXt( int lay, int iEntr, int lr, TGraph* gr ) {
372 unsigned vsize = m_vt.size();
373 if ( vsize < 15 ) return false;
374
375 double tmax = m_vt[vsize - 1];
376 double tm0 = 300.0; // 280->300
377 double tm1 = 500.0;
378 double tm2 = 700.0;
379 if ( lay < 8 )
380 {
381 tm0 = 200.0; // 180->200
382 tm1 = 300.0;
383 tm2 = 540.0;
384 }
385 int n0 = 0;
386 for ( unsigned i = 0; i < vsize; i++ )
387 {
388 if ( m_vt[i] < tm0 ) n0++;
389 }
390
391 int nCut = 30; // 25
392 if ( lay < 8 ) nCut = 20; // 15
393
394 double entries1 = 0.0;
395 double entries2 = 0.0;
396 for ( unsigned i = 0; i < vsize; i++ )
397 {
398 entries1 += m_entries[i];
399 if ( m_vt[i] < 150. ) entries2 += m_entries[i];
400 }
401 if ( ( entries1 * 0.9 ) < entries2 ) return false;
402
403 if ( n0 < nCut ) return false;
404 if ( tmax < tm1 ) return funXt0( lay, iEntr, lr, gr );
405 else return funXt1( lay, iEntr, lr, gr );
406}
407
408bool XtInteCalib::funXt0( int lay, int iEntr, int lr, TGraph* gr ) {
409 double tCut = 300.0;
410 if ( lay < 8 ) tCut = 200.0;
411
412 int npoint = 0;
413 if ( m_vt[0] > 0.0 )
414 {
415 gr->SetPoint( 0, 0, 0 );
416 npoint++;
417 }
418 for ( unsigned i = 0; i < m_vt.size(); i++ )
419 {
420 double delt = 0.0;
421 if ( i > 10 ) delt = m_vt[i] - m_vt[i - 1];
422 if ( ( i > 10 ) && ( ( delt > 100. ) || ( ( delt > 60. ) && ( m_vt[i - 1] < tCut ) ) ) )
423 break;
424
425 if ( m_vt[i] < 0.0 ) gr->SetPoint( npoint, m_vt[i], 0.0 );
426 else gr->SetPoint( npoint, m_vt[i], m_vd[i] );
427 npoint++;
428 }
429 return true;
430}
431
432bool XtInteCalib::funXt1( int lay, int iEntr, int lr, TGraph* gr ) {
433 unsigned vsize = m_vt.size();
434 double tmax = m_vt[vsize - 1];
435 double t1;
436 if ( lay < 8 )
437 {
438 if ( tmax < 540 ) t1 = 300.;
439 else t1 = 540.;
440 }
441 else
442 {
443 if ( tmax < 700 ) t1 = 500.;
444 else t1 = 660.;
445 }
446
447 bool fgfit = false;
448 int np = 0;
449 double c0 = 0.0;
450 double c1 = 0.0;
451 TGraph* grPol1 = new TGraph();
452 for ( unsigned i = 0; i < vsize; i++ )
453 {
454 if ( m_vt[i] >= t1 )
455 {
456 grPol1->SetPoint( np, m_vt[i], m_vd[i] );
457 np++;
458 }
459 }
460 double t2;
461 if ( np < 5 ) { t2 = t1; }
462 else
463 {
464 double x2;
465 grPol1->GetPoint( 0, t2, x2 );
466 grPol1->Fit( "pol1", "Q", "", t2, m_vt[vsize - 1] );
467 c0 = grPol1->GetFunction( "pol1" )->GetParameter( 0 );
468 c1 = grPol1->GetFunction( "pol1" )->GetParameter( 1 );
469
470 if ( c1 < 0 )
471 {
472 grPol1->Fit( "pol0", "Q", "", t2, m_vt[vsize - 1] );
473 c0 = grPol1->GetFunction( "pol0" )->GetParameter( 0 );
474 c1 = 0.0;
475 }
476 fgfit = true;
477 }
478
479 double tCut = 300.0;
480 if ( lay < 8 ) tCut = 200.0;
481
482 int npoint = 0;
483 if ( m_vt[0] > 0.0 )
484 {
485 gr->SetPoint( 0, 0, 0 );
486 npoint++;
487 }
488 for ( unsigned i = 0; i < vsize; i++ )
489 {
490 double delt = 0.0;
491 if ( i > 10 ) delt = m_vt[i] - m_vt[i - 1];
492 if ( ( i > 10 ) && ( ( delt > 100. ) || ( ( delt > 60. ) && ( m_vt[i - 1] < tCut ) ) ) )
493 break;
494
495 if ( m_vt[i] < t2 )
496 {
497 if ( m_vt[i] < 0.0 ) gr->SetPoint( npoint, m_vt[i], 0.0 );
498 else gr->SetPoint( npoint, m_vt[i], m_vd[i] );
499 npoint++;
500 }
501 else
502 {
503 double dist;
504 if ( !fgfit ) dist = m_vd[i];
505 else dist = c1 * m_vt[i] + c0;
506 gr->SetPoint( npoint, m_vt[i], dist );
507 npoint++;
508 }
509 }
510 return true;
511}
512
513int XtInteCalib::findXtEntr( int lay, int iEntr, int lr ) const {
514 int id0 = 8;
515 int id1 = 9;
516 int idmax = 17;
517 int entrId = -1;
518 if ( iEntr <= id0 )
519 {
520 int id = -1;
521 for ( int i = iEntr; i <= id0; i++ )
522 {
523 if ( m_fgFit[lay][i][lr] )
524 {
525 id = i;
526 break;
527 }
528 }
529 if ( -1 != id ) entrId = id;
530 else
531 {
532 for ( int i = iEntr; i >= 0; i-- )
533 {
534 if ( m_fgFit[lay][i][lr] )
535 {
536 id = i;
537 break;
538 }
539 }
540 if ( -1 != id ) entrId = id;
541 else
542 {
543 for ( int i = id1; i <= idmax; i++ )
544 {
545 if ( m_fgFit[lay][i][lr] )
546 {
547 id = i;
548 break;
549 }
550 }
551 entrId = id;
552 }
553 }
554 }
555 else
556 {
557 int id = -1;
558 for ( int i = iEntr; i >= id1; i-- )
559 {
560 if ( m_fgFit[lay][i][lr] )
561 {
562 id = i;
563 break;
564 }
565 }
566 if ( -1 != id ) entrId = id;
567 else
568 {
569 for ( int i = iEntr; i < idmax; i++ )
570 {
571 if ( m_fgFit[lay][i][lr] )
572 {
573 id = i;
574 break;
575 }
576 }
577 if ( -1 != id ) entrId = id;
578 else
579 {
580 for ( int i = id1; i >= 0; i-- )
581 {
582 if ( m_fgFit[lay][i][lr] )
583 {
584 id = i;
585 break;
586 }
587 }
588 entrId = id;
589 }
590 }
591 }
592 if ( -1 == entrId )
593 {
594 cout << "find EntrId error "
595 << "layer " << lay << " iEntr " << iEntr << " lr " << lr << endl;
596 }
597
598 return entrId;
599}
600
601int XtInteCalib::findXtEntrEdge( int lay, int iEntr, int lr ) const {
602 int id0 = 8;
603 int id1 = 9;
604 int idmax = 17;
605 int entrId = -1;
606 if ( iEntr <= id0 )
607 {
608 int id = -1;
609 for ( int i = iEntr; i <= id0; i++ )
610 {
611 if ( m_fgEdge[lay][i][lr] )
612 {
613 id = i;
614 break;
615 }
616 }
617 if ( -1 != id ) entrId = id;
618 else
619 {
620 for ( int i = iEntr; i >= 0; i-- )
621 {
622 if ( m_fgEdge[lay][i][lr] )
623 {
624 id = i;
625 break;
626 }
627 }
628 if ( -1 != id ) entrId = id;
629 else
630 {
631 for ( int i = id1; i <= idmax; i++ )
632 {
633 if ( m_fgEdge[lay][i][lr] )
634 {
635 id = i;
636 break;
637 }
638 }
639 entrId = id;
640 }
641 }
642 }
643 else
644 {
645 int id = -1;
646 for ( int i = iEntr; i >= id1; i-- )
647 {
648 if ( m_fgEdge[lay][i][lr] )
649 {
650 id = i;
651 break;
652 }
653 }
654 if ( -1 != id ) entrId = id;
655 else
656 {
657 for ( int i = iEntr; i < idmax; i++ )
658 {
659 if ( m_fgEdge[lay][i][lr] )
660 {
661 id = i;
662 break;
663 }
664 }
665 if ( -1 != id ) entrId = id;
666 else
667 {
668 for ( int i = id1; i >= 0; i-- )
669 {
670 if ( m_fgEdge[lay][i][lr] )
671 {
672 id = i;
673 break;
674 }
675 }
676 entrId = id;
677 }
678 }
679 }
680 if ( -1 == entrId )
681 {
682 cout << "find EntrId error for cell edge "
683 << "layer " << lay << " iEntr " << iEntr << " lr " << lr << endl;
684 }
685
686 return entrId;
687}
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)
TGraph * gr
#define dmax(a, b)
DOUBLE_PRECISION tr[3]
int n2
Definition SD0Tag.cxx:59
int n1
Definition SD0Tag.cxx:58
virtual void init(TObjArray *hlist, MdcCosGeom *pGeom)=0
Definition CalibBase.cpp:12
virtual void mergeHist(TFile *fhist)=0
virtual void calib(MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)=0
void calib(MdcCalibConst *calconst, TObjArray *newXtList, TObjArray *r2tList)
void mergeHist(TFile *fhist)
void init(TObjArray *hlist, MdcCosGeom *pGeom)
int t()
Definition t.c:1