BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
IniMdcCalib.cxx
Go to the documentation of this file.
2
3#include "GaudiKernel/Bootstrap.h"
4#include "GaudiKernel/IMessageSvc.h"
5#include "GaudiKernel/ISvcLocator.h"
6#include "GaudiKernel/MsgStream.h"
7#include "GaudiKernel/StatusCode.h"
8
9#include "GaudiKernel/IDataProviderSvc.h"
10#include "GaudiKernel/PropertyMgr.h"
11#include "GaudiKernel/SmartDataPtr.h"
12
13#include "EvTimeEvent/RecEsTime.h"
14#include "EventModel/Event.h"
15#include "EventModel/EventHeader.h"
16#include "Identifier/Identifier.h"
17#include "Identifier/MdcID.h"
18#include "MdcRawEvent/MdcDigi.h"
19#include "TStyle.h"
20
21#include <cmath>
22#include <cstring>
23#include <fstream>
24#include <iomanip>
25#include <iostream>
26#include <string>
27
28#include "TF1.h"
29
30using namespace Event;
31
33
35
37 int iEs;
38 for ( int lay = 0; lay < MdcCalNLayer; lay++ )
39 {
40 delete m_hlaymapT[lay];
41 delete m_htdc[lay];
42 delete m_htraw[lay];
43 for ( iEs = 0; iEs < m_param.nEsFlag; iEs++ )
44 {
45 delete m_htdcTes[lay][iEs];
46 delete m_htrawTes[lay][iEs];
47 }
48
49 delete m_hlaymapQ[lay];
50 delete m_hqraw[lay];
51 }
52
53 for ( int wir = 0; wir < MdcCalTotCell; wir++ )
54 {
55 delete m_htrawCel[wir];
56 delete m_hqrawCel[wir];
57 }
58
59 delete m_hLayerHitmapT;
60 delete m_hWireHitMapT;
61
62 delete m_hLayerHitmapQ;
63 delete m_hWireHitMapQ;
64 for ( iEs = 0; iEs < m_param.nEsFlag; iEs++ ) delete m_hTes[iEs];
65 delete m_hTesAllFlag;
66 delete m_hTesAll;
67 delete m_hTesCal;
68 delete m_hTesFlag;
69
70 delete m_fdcom;
71 delete m_fdTmap;
72 delete m_fdTraw;
73 delete m_fdTrawCel;
74 delete m_fdTrawTes;
75 delete m_fdQmap;
76 delete m_fdQraw;
77 delete m_fdQrawCel;
78}
79
80void IniMdcCalib::initialize( TObjArray* hlist, IMdcGeomSvc* mdcGeomSvc,
81 IMdcCalibFunSvc* mdcFunSvc, IMdcUtilitySvc* mdcUtilitySvc ) {
82 IMessageSvc* msgSvc;
83 Gaudi::svcLocator()->service( "MessageSvc", msgSvc );
84 MsgStream log( msgSvc, "IniMdcCalib" );
85 log << MSG::INFO << "IniMdcCalib::initialize()" << endmsg;
86
87 m_hlist = hlist;
88 m_mdcGeomSvc = mdcGeomSvc;
89 m_mdcFunSvc = mdcFunSvc;
90 m_mdcUtilitySvc = mdcUtilitySvc;
91
92 int lay;
93 int cel;
94 int wir;
95 int ncel;
96 char hname[200];
97
98 m_nWire = m_mdcGeomSvc->getWireSize();
99 m_nLayer = m_mdcGeomSvc->getLayerSize();
100
101 m_fdcom = new TFolder( "Common", "Common" );
102 m_hlist->Add( m_fdcom );
103
104 m_fdTmap = new TFolder( "Thitmap", "Thitmap" );
105 m_hlist->Add( m_fdTmap );
106
107 m_fdTraw = new TFolder( "Traw", "Traw" );
108 m_hlist->Add( m_fdTraw );
109
110 m_fdTrawCel = new TFolder( "TrawCell", "TrawCell" );
111 m_hlist->Add( m_fdTrawCel );
112
113 m_fdTrawTes = new TFolder( "TrawTes", "TrawTes" );
114 m_hlist->Add( m_fdTrawTes );
115
116 m_fdQmap = new TFolder( "Qhitmap", "Qhitmap" );
117 m_hlist->Add( m_fdQmap );
118
119 m_fdQraw = new TFolder( "Qraw", "Qraw" );
120 m_hlist->Add( m_fdQraw );
121
122 m_fdQrawCel = new TFolder( "QrawCell", "QrawCell" );
123 m_hlist->Add( m_fdQrawCel );
124
125 m_hLayerHitmapT = new TH1F( "T_Hitmap_Layer", "", 43, -0.5, 42.5 );
126 m_fdcom->Add( m_hLayerHitmapT );
127
128 m_hWireHitMapT = new TH1F( "T_Hitmap_Wire", "", 6796, -0.5, 6795.5 );
129 m_fdcom->Add( m_hWireHitMapT );
130
131 m_hLayerHitmapQ = new TH1F( "Q_Hitmap_Layer", "", 43, -0.5, 42.5 );
132 m_fdcom->Add( m_hLayerHitmapQ );
133
134 m_hWireHitMapQ = new TH1F( "Q_Hitmap_Wire", "", 6796, -0.5, 6795.5 );
135 m_fdcom->Add( m_hWireHitMapQ );
136
137 int iEs;
138 for ( iEs = 0; iEs < m_param.nEsFlag; iEs++ )
139 {
140 sprintf( hname, "Tes_%d", m_param.esFlag[iEs] );
141 m_hTes[iEs] = new TH1F( hname, "", 2000, 0, 2000 );
142 m_fdcom->Add( m_hTes[iEs] );
143 }
144
145 m_hTesAllFlag = new TH1F( "TesAllFlag", "", 300, -0.5, 299.5 );
146 m_fdcom->Add( m_hTesAllFlag );
147
148 m_hTesAll = new TH1F( "TesAll", "", 2000, 0, 2000 );
149 m_fdcom->Add( m_hTesAll );
150
151 m_hTesCal = new TH1F( "TesCal", "", 2000, 0, 2000 );
152 m_fdcom->Add( m_hTesCal );
153
154 m_hTesFlag = new TH1F( "Tes_Flag", "", 300, -0.5, 299.5 );
155 m_fdcom->Add( m_hTesFlag );
156
157 int tdcNbin = 1000;
158 double tdcmin = 0;
159 double tdcmax = 20000;
160
161 int trawNbin = 1000;
162 double trawmin = 0;
163 double trawmax = 2000;
164
165 for ( lay = 0; lay < m_nLayer; lay++ )
166 {
167 ncel = m_mdcGeomSvc->Layer( lay )->NCell();
168
169 sprintf( hname, "T_hitmap_Lay%02d", lay );
170 m_hlaymapT[lay] = new TH1F( hname, "", ncel, -0.5, (float)ncel - 0.5 );
171 m_fdTmap->Add( m_hlaymapT[lay] );
172
173 sprintf( hname, "TDC_Lay%02d", lay );
174 m_htdc[lay] = new TH1F( hname, "", tdcNbin, tdcmin, tdcmax );
175 m_fdTraw->Add( m_htdc[lay] );
176
177 sprintf( hname, "Traw_Lay%02d", lay );
178 m_htraw[lay] = new TH1F( hname, "", trawNbin, trawmin, trawmax );
179 m_fdTraw->Add( m_htraw[lay] );
180
181 for ( iEs = 0; iEs < m_param.nEsFlag; iEs++ )
182 {
183 sprintf( hname, "TDC_Lay%02d_Tes%d", lay, m_param.esFlag[iEs] );
184 m_htdcTes[lay][iEs] = new TH1F( hname, "", tdcNbin, tdcmin, tdcmax );
185 m_fdTrawTes->Add( m_htdcTes[lay][iEs] );
186
187 sprintf( hname, "Traw_Lay%02d_Tes%d", lay, m_param.esFlag[iEs] );
188 m_htrawTes[lay][iEs] = new TH1F( hname, "", trawNbin, trawmin, trawmax );
189 m_fdTrawTes->Add( m_htrawTes[lay][iEs] );
190 }
191
192 sprintf( hname, "Q_hitmap_Lay%02d", lay );
193 m_hlaymapQ[lay] = new TH1F( hname, "", ncel, -0.5, (float)ncel - 0.5 );
194 m_fdQmap->Add( m_hlaymapQ[lay] );
195
196 sprintf( hname, "Qraw_Lay%02d", lay );
197 m_hqraw[lay] = new TH1F( hname, "", 2000, 0, 4000 );
198 m_fdQraw->Add( m_hqraw[lay] );
199 }
200
201 for ( wir = 0; wir < MdcCalTotCell; wir++ )
202 {
203 lay = m_mdcGeomSvc->Wire( wir )->Layer();
204 cel = m_mdcGeomSvc->Wire( wir )->Cell();
205
206 sprintf( hname, "Traw_%02d_%03d_%04d", lay, cel, wir );
207 m_htrawCel[wir] = new TH1F( hname, "", trawNbin, trawmin, trawmax );
208 m_fdTrawCel->Add( m_htrawCel[wir] );
209
210 sprintf( hname, "Qraw_%02d_%03d_%04d", lay, cel, wir );
211 m_hqrawCel[wir] = new TH1F( hname, "", 2000, 0, 4000 );
212 m_fdQrawCel->Add( m_hqrawCel[wir] );
213 }
214}
215
217 IMessageSvc* msgSvc;
218 Gaudi::svcLocator()->service( "MessageSvc", msgSvc );
219 MsgStream log( msgSvc, "IniMdcCalib" );
220 log << MSG::DEBUG << "IniMdcCalib::fillHist()" << endmsg;
221
222 int lay;
223 int cel;
224 int wir;
225 int digiId;
226 unsigned fgOverFlow;
227 double tdc;
228 double traw;
229 double adc;
230 double qraw;
231 Identifier id;
232
233 IDataProviderSvc* eventSvc = NULL;
234 Gaudi::svcLocator()->service( "EventDataSvc", eventSvc );
235
236 SmartDataPtr<Event::EventHeader> eventHeader( eventSvc, "/Event/EventHeader" );
237 if ( !eventHeader )
238 {
239 log << MSG::FATAL << "Could not find Event Header" << endmsg;
240 return ( StatusCode::FAILURE );
241 }
242 int iRun = eventHeader->runNumber();
243 int iEvt = eventHeader->eventNumber();
244
245 // get EsTimeCol
246 double tes = -9999.0;
247 int esTimeflag = -1;
248 SmartDataPtr<RecEsTimeCol> aevtimeCol( eventSvc, "/Event/Recon/RecEsTimeCol" );
249 if ( ( !aevtimeCol ) || ( aevtimeCol->size() == 0 ) )
250 {
251 tes = -9999.0;
252 esTimeflag = -1;
253 }
254 else
255 {
256 RecEsTimeCol::iterator iter_evt = aevtimeCol->begin();
257 for ( ; iter_evt != aevtimeCol->end(); iter_evt++ )
258 {
259 tes = ( *iter_evt )->getTest();
260 esTimeflag = ( *iter_evt )->getStat();
261 // cout << "tes " << tes << endl;
262 }
263 }
264 m_hTesAllFlag->Fill( esTimeflag );
265 m_hTesAll->Fill( tes );
266 m_hTesFlag->Fill( esTimeflag );
267 int nTesFlag = -1;
268 for ( int iEs = 0; iEs < m_param.nEsFlag; iEs++ )
269 {
270 if ( esTimeflag == m_param.esFlag[iEs] )
271 {
272 m_hTes[iEs]->Fill( tes );
273 nTesFlag = iEs;
274 break;
275 }
276 }
277 bool fgFillTes = false;
278 if ( !m_param.esCut ) { fgFillTes = true; }
279 else if ( ( nTesFlag > -1 ) && ( tes > m_param.tesMin ) && ( tes < m_param.tesMax ) )
280 { fgFillTes = true; }
281
282 if ( fgFillTes ) m_hTesCal->Fill( tes );
283
284 // cout << setw(10) << iRun << setw(10) << iEvt << setw(10) << tes << endl;
285
286 // retrieve Mdc digi
287 SmartDataPtr<MdcDigiCol> mdcDigiCol( eventSvc, "/Event/Digi/MdcDigiCol" );
288 if ( !mdcDigiCol ) { log << MSG::FATAL << "Could not find event" << endmsg; }
289
290 MdcDigiCol::iterator iter = mdcDigiCol->begin();
291 digiId = 0;
292 for ( ; iter != mdcDigiCol->end(); iter++, digiId++ )
293 {
294 MdcDigi* aDigi = ( *iter );
295 id = ( aDigi )->identify();
296
297 lay = MdcID::layer( id );
298 cel = MdcID::wire( id );
299 wir = m_mdcGeomSvc->Wire( lay, cel )->Id();
300
301 tdc = ( aDigi )->getTimeChannel();
302 adc = ( aDigi )->getChargeChannel();
303 fgOverFlow = ( aDigi )->getOverflow();
304
305 if ( ( ( fgOverFlow & 1 ) != 0 ) || ( ( fgOverFlow & 12 ) != 0 ) ||
306 ( aDigi->getTimeChannel() == 0x7FFFFFFF ) ||
307 ( aDigi->getChargeChannel() == 0x7FFFFFFF ) )
308 continue;
309
310 traw = tdc * MdcCalTdcCnv;
311 qraw = adc * MdcCalAdcCnv;
312
313 m_hLayerHitmapT->Fill( lay );
314 m_hWireHitMapT->Fill( wir );
315 m_hlaymapT[lay]->Fill( cel );
316
317 m_hLayerHitmapQ->Fill( lay );
318 m_hWireHitMapQ->Fill( wir );
319 m_hlaymapQ[lay]->Fill( cel );
320
321 if ( fgFillTes )
322 {
323 traw -= tes;
324 traw += m_param.timeShift;
325
326 m_htdc[lay]->Fill( tdc );
327 m_htraw[lay]->Fill( traw );
328 m_hqraw[lay]->Fill( qraw );
329
330 if ( nTesFlag > -1 )
331 {
332 m_htdcTes[lay][nTesFlag]->Fill( tdc );
333 m_htrawTes[lay][nTesFlag]->Fill( traw );
334 }
335
336 m_htrawCel[wir]->Fill( traw );
337 m_hqrawCel[wir]->Fill( qraw );
338 }
339 }
340 return 1;
341}
342
344
346 IMessageSvc* msgSvc;
347 Gaudi::svcLocator()->service( "MessageSvc", msgSvc );
348 MsgStream log( msgSvc, "IniMdcCalib" );
349 log << MSG::DEBUG << "IniMdcCalib::updateConst()" << endmsg;
350
351 int lay;
352 int wir;
353 double t0Fit[MdcCalNLayer];
354 double t0Cal[MdcCalNLayer];
355 double tmax[MdcCalNLayer];
356 double initT0 = m_param.initT0;
357
358 int fitTminFg[MdcCalNLayer];
359 int fitTmaxFg[MdcCalNLayer];
360 double chisq;
361 double ndf;
362 double chindfTmin[MdcCalNLayer];
363 double chindfTmax[MdcCalNLayer];
364 char funname[200];
365
366 // fit Tmin
367 TF1* ftmin[MdcCalNLayer];
368 // sprintf(funname, "ftmin", lay);
369 // TF1* ftmin = new TF1(funname, funTmin, 0, 150, 6);
370 for ( lay = 0; lay < m_nLayer; lay++ )
371 {
372 fitTminFg[lay] = 0;
373 chindfTmin[lay] = -1;
374 sprintf( funname, "ftmin%02d", lay );
375 ftmin[lay] = new TF1( funname, funTmin, 0, 150, 6 );
376
377 if ( 1 == m_param.fgCalib[lay] )
378 {
379 Stat_t nEntryTot = 0;
380 for ( int ibin = 1; ibin <= 25; ibin++ )
381 {
382 Stat_t entry = m_htraw[lay]->GetBinContent( ibin );
383 nEntryTot += entry;
384 }
385 double c0Ini = (double)nEntryTot / 25.0;
386 double c1Ini = ( m_htraw[lay]->GetMaximum() ) - c0Ini;
387
388 ftmin[lay]->SetParameter( 0, c0Ini );
389 ftmin[lay]->SetParameter( 1, c1Ini );
390 ftmin[lay]->SetParameter( 2, 0 );
391 ftmin[lay]->SetParameter( 4, initT0 );
392 ftmin[lay]->SetParameter( 5, 2 );
393
394 m_htraw[lay]->Fit( funname, "Q", "", m_param.tminFitRange[lay][0],
395 m_param.tminFitRange[lay][1] );
396 chisq = ftmin[lay]->GetChisquare();
397 gStyle->SetOptFit( 11 );
398 ndf = ftmin[lay]->GetNDF();
399 chindfTmin[lay] = chisq / ndf;
400 if ( chindfTmin[lay] < m_param.tminFitChindf )
401 {
402 fitTminFg[lay] = 1;
403 t0Fit[lay] = ftmin[lay]->GetParameter( 4 );
404 t0Fit[lay] += m_param.t0Shift;
405 t0Cal[lay] = t0Fit[lay] - m_param.timeShift;
406 }
407 }
408
409 if ( 0 == fitTminFg[lay] )
410 {
411 wir = m_mdcGeomSvc->Wire( lay, 0 )->Id();
412 t0Cal[lay] = calconst->getT0( wir );
413 t0Fit[lay] = t0Cal[lay] + m_param.timeShift;
414 }
415 }
416
417 // fit Tmax
418 TF1* ftmax[MdcCalNLayer];
419 for ( lay = 0; lay < m_nLayer; lay++ )
420 {
421 fitTmaxFg[lay] = 0;
422 chindfTmax[lay] = -1;
423 sprintf( funname, "ftmax%02d", lay );
424 ftmax[lay] = new TF1( funname, funTmax, 250, 500, 4 );
425
426 if ( 1 == m_param.fgCalib[lay] )
427 {
428 ftmax[lay]->SetParameter( 2, m_param.initTm[lay] );
429 ftmax[lay]->SetParameter( 3, 10 );
430
431 m_htraw[lay]->Fit( funname, "Q+", "", m_param.tmaxFitRange[lay][0],
432 m_param.tmaxFitRange[lay][1] );
433 gStyle->SetOptFit( 11 );
434 chisq = ftmax[lay]->GetChisquare();
435 ndf = ftmax[lay]->GetNDF();
436 chindfTmax[lay] = chisq / ndf;
437 if ( chindfTmax[lay] < m_param.tmaxFitChindf )
438 {
439 fitTmaxFg[lay] = 1;
440 tmax[lay] = ftmax[lay]->GetParameter( 2 );
441 }
442 }
443
444 if ( 0 == fitTmaxFg[lay] )
445 { tmax[lay] = ( calconst->getXtpar( lay, 0, 0, 6 ) ) + t0Fit[lay]; }
446 }
447
448 // output for check
449 ofstream ft0( "iniT0.dat" );
450 for ( lay = 0; lay < m_nLayer; lay++ )
451 {
452 ft0 << setw( 5 ) << lay << setw( 3 ) << fitTminFg[lay] << setw( 12 ) << t0Cal[lay]
453 << setw( 12 ) << t0Fit[lay] << setw( 12 ) << chindfTmin[lay] << setw( 5 )
454 << fitTmaxFg[lay] << setw( 12 ) << tmax[lay] << setw( 12 ) << tmax[lay] - t0Fit[lay]
455 << setw( 12 ) << chindfTmax[lay] << endl;
456 }
457 ft0.close();
458 cout << "iniT0.dat was written." << endl;
459
460 // set T0
461 int i;
462 int nwire = m_mdcGeomSvc->getWireSize();
463 for ( i = 0; i < nwire; i++ )
464 {
465 lay = m_mdcGeomSvc->Wire( i )->Layer();
466 if ( 1 == m_param.fgCalib[lay] )
467 {
468 calconst->resetT0( i, t0Cal[lay] );
469 calconst->resetDelT0( i, 0.0 );
470 }
471 }
472
473 if ( 0 == m_param.fgIniCalConst )
474 {
475 // set X-T
476 int lr;
477 int iEntr;
478 int ord;
479 double xtpar;
480 double xtini[8] = { 0, 0.03, 0, 0, 0, 0, 999.9, 0 };
481 for ( lay = 0; lay < m_nLayer; lay++ )
482 {
483 if ( 1 != m_param.fgCalib[lay] ) continue;
484
485 for ( iEntr = 0; iEntr < MdcCalNENTRXT; iEntr++ )
486 {
487 for ( lr = 0; lr < MdcCalLR; lr++ )
488 {
489 for ( ord = 0; ord < MdcCalXtNPars; ord++ )
490 {
491 if ( 6 == ord ) { xtpar = tmax[lay] - t0Fit[lay]; }
492 else { xtpar = xtini[ord]; }
493 calconst->resetXtpar( lay, iEntr, lr, ord, xtpar );
494 }
495 }
496 }
497 }
498
499 // set Q-T
500 for ( lay = 0; lay < m_nLayer; lay++ )
501 {
502 calconst->resetQtpar0( lay, 0.0 );
503 calconst->resetQtpar1( lay, 0.0 );
504 }
505
506 // set S-D
507 int bin;
508 double sdpar = m_param.initSigma; // mm
509 for ( lay = 0; lay < m_nLayer; lay++ )
510 {
511 for ( iEntr = 0; iEntr < MdcCalNENTRSD; iEntr++ )
512 {
513 for ( lr = 0; lr < 2; lr++ )
514 {
515 for ( bin = 0; bin < MdcCalSdNBIN; bin++ )
516 { calconst->resetSdpar( lay, iEntr, lr, bin, sdpar ); }
517 }
518 }
519 }
520 }
521 else if ( 2 == m_param.fgIniCalConst )
522 {
523 int lr;
524 int iEntr;
525 double xtpar;
526 for ( lay = 0; lay < m_nLayer; lay++ )
527 {
528 for ( iEntr = 0; iEntr < MdcCalNENTRXT; iEntr++ )
529 {
530 for ( lr = 0; lr < MdcCalLR; lr++ )
531 {
532 xtpar = tmax[lay] - t0Fit[lay];
533 calconst->resetXtpar( lay, iEntr, lr, 6, xtpar );
534 }
535 }
536 }
537 }
538
539 // delete ftmin;
540 for ( lay = 0; lay < MdcCalNLayer; lay++ )
541 {
542 delete ftmin[lay];
543 delete ftmax[lay];
544 }
545 return 1;
546}
547
548Double_t IniMdcCalib::funTmin( Double_t* x, Double_t* par ) {
549 Double_t fitval;
550 fitval = par[0] + par[1] * exp( -par[2] * ( x[0] - par[3] ) ) /
551 ( 1 + exp( -( x[0] - par[4] ) / par[5] ) );
552 return fitval;
553}
554
555Double_t IniMdcCalib::funTmax( Double_t* x, Double_t* par ) {
556 Double_t fitval;
557 fitval = par[0] + par[1] / ( 1 + exp( ( x[0] - par[2] ) / par[3] ) );
558 return fitval;
559}
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)
Double_t x[10]
EvtComplex exp(const EvtComplex &c)
EvtStreamInputIterator< typename Generator::result_type > iter(Generator gen, int N=0)
*******INTEGER m_nBinMax INTEGER m_NdiMax !No of bins in histogram for cell exploration division $ !Last vertex $ !Last active cell $ !Last cell in buffer $ !No of sampling when dividing cell $ !No of function total $ !Flag for random ceel for $ !Flag for type of for WtMax $ !Flag which decides whether vertices are included in the sampling $ entire domain is hyp !Maximum effective eevents per bin
Definition FoamA.h:85
const int MdcCalNLayer
Definition MdcCalParams.h:6
const int MdcCalNENTRSD
const double MdcCalTdcCnv
const double MdcCalAdcCnv
const int MdcCalXtNPars
const int MdcCalSdNBIN
const int MdcCalNENTRXT
const int MdcCalTotCell
Definition MdcCalParams.h:9
const int MdcCalLR
IMessageSvc * msgSvc()
int updateConst(MdcCalibConst *calconst)
void printCut() const
int fillHist(MdcCalEvent *event)
void initialize(TObjArray *hlist, IMdcGeomSvc *mdcGeomSvc, IMdcCalibFunSvc *mdcFunSvc, IMdcUtilitySvc *mdcUtilitySvc)
void resetQtpar0(int lay, double val)
void resetSdpar(int lay, int entr, int lr, int bin, double val)
void resetXtpar(int lay, int entr, int lr, int order, double val)
void resetQtpar1(int lay, double val)
double getT0(int wireid) const
void resetDelT0(int wireid, double val)
double getXtpar(int lay, int entr, int lr, int order)
void resetT0(int wireid, double val)
static int layer(const Identifier &id)
Values of different levels (failure returns 0).
Definition MdcID.cxx:47
static int wire(const Identifier &id)
Definition MdcID.cxx:52
unsigned int getChargeChannel() const
Definition RawData.cxx:35
unsigned int getTimeChannel() const
Definition RawData.cxx:32