StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StSvtBadAnodesMaker.cxx
1 /***************************************************************************
2  *
3  * $Id: StSvtBadAnodesMaker.cxx,v 1.2 2004/02/04 16:13:56 munhoz Exp $
4  *
5  * Author: Marcelo Munhoz
6  ***************************************************************************
7  *
8  * Description: SVT Bad Anodes Maker
9  *
10  ***************************************************************************
11  *
12  * $Log: StSvtBadAnodesMaker.cxx,v $
13  * Revision 1.2 2004/02/04 16:13:56 munhoz
14  * fix few problems with pedestal reading
15  *
16  * Revision 1.1 2004/02/03 20:06:28 munhoz
17  * first version of bad anode maker
18  *
19  *
20  **************************************************************************/
21 
22 #include <fstream>
23 
24 #include "TH2.h"
25 #include "TString.h"
26 #include "TFile.h"
27 #include "TKey.h"
28 #include "TObjectSet.h"
29 #include "StMessMgr.h"
30 
31 #include "StSvtBadAnodesMaker.h"
32 #include "StarClassLibrary/StSequence.hh"
33 #include "StSvtClassLibrary/StSvtHybridBadAnodes.hh"
34 #include "StSvtClassLibrary/StSvtHybridData.hh"
35 #include "StSvtClassLibrary/StSvtHybridCollection.hh"
36 #include "StSvtClassLibrary/StSvtData.hh"
37 #include "StSvtClassLibrary/StSvtHybridPed.hh"
38 #include "map.h"
39 
40 #define MAX_NUMBER_OF_ANODES 240
41 
42 ClassImp(StSvtBadAnodesMaker)
43 
44 //_____________________________________________________________________________
45 StSvtBadAnodesMaker::StSvtBadAnodesMaker(const char *name):StMaker(name)
46 {
47  mSvtData = NULL;
48  mHybridData = NULL;
49  mSvtBadAnodes = NULL;
50  mHybridBadAnodes = NULL;
51  mSvtPed = NULL;
52  mSvtRMSPed = NULL;
53 
54  mBadAnodesSet = NULL;
55 
56  mRmsScaleFactor = 16;
57 
58  NULL_ADC = 0; //
59  OVERLOADED_ADC = 255; //
60  BAD_RMS = 10; // in ADC counts
61  BAD_MEAN_PED_MIN = 0; //
62  BAD_MEAN_PED_MAX = 255; //
63  BAD_MEAN_RMS_MIN = 0; //
64  BAD_MEAN_RMS_MAX = 255; //
65 
66  NULL_ADC_THRESHOLD = 129; //
67  OVERLOADED_ADC_THRESHOLD = 129; // in number of time bins
68  OCCUP_THRESHOLD = 129; //
69  RMS_THRESHOLD = 129; //
70 
71  FREQ_OVERLOADED_ADC = 1.1; //
72  FREQ_NULL_ADC = 1.1; // fraction of events
73  FREQ_OCCUP = 1.1; //
74 
75  mFileName = new TString("badAnodes.root");
76 }
77 
78 //_____________________________________________________________________________
79 StSvtBadAnodesMaker::~StSvtBadAnodesMaker()
80 {}
81 
82 void StSvtBadAnodesMaker::setOutputFile(const char* name)
83 {
84  if (mFileName)
85  delete mFileName;
86 
87  mFileName = new TString(name);
88 }
89 
90 //_____________________________________________________________________________
91 Int_t StSvtBadAnodesMaker::Init()
92 {
93  if (Debug()) gMessMgr->Debug() << "StSvtBadAnodesMaker::Init" << endm;
94 
95  mEvents = 0;
96 
97  setSvtData();
98  setPedestal();
99  setRMSPedestal();
100 
101  setSvtBadAnodes();
102 
103  /*
104  setPedestal();
105 
106  setRMSPedestal();
107  */
108 
109  if (Debug()) bookHistograms();
110 
111  gMessMgr->Info() << "NULL_ADC = " << NULL_ADC << endm;
112  gMessMgr->Info() << "OVERLOADED_ADC = " << OVERLOADED_ADC << endm;
113 
114  gMessMgr->Info() << "THRESHOLD NULL_ADC (fraction of time bins) = " << NULL_ADC_THRESHOLD << endm;
115  gMessMgr->Info() << "THRESHOLD OVERLOADED_ADC (fraction of time bins) = " << OVERLOADED_ADC_THRESHOLD << endm;
116  gMessMgr->Info() << "THRESHOLD OCCUP (fraction of time bins) = " << OCCUP_THRESHOLD << endm;
117 
118  gMessMgr->Info() << "FREQUENCY OVERLOADED_ADC (fraction of events) = " << FREQ_OVERLOADED_ADC << endm;
119  gMessMgr->Info() << "FREQUENCY NULL_ADC (fraction of events) = " << FREQ_NULL_ADC << endm;
120  gMessMgr->Info() << "FREQUENCY OCCUP (fraction of events) = " << FREQ_OCCUP << endm;
121 
122  gMessMgr->Info() << "BAD MEAN PED_MIN = " << BAD_MEAN_PED_MIN << endm;
123  gMessMgr->Info() << "BAD MEAN PED_MAX = " << BAD_MEAN_PED_MAX << endm;
124  gMessMgr->Info() << "BAD MEAN RMS_MIN = " << BAD_MEAN_RMS_MIN << endm;
125  gMessMgr->Info() << "BAD MEAN RMS_MAX = " << BAD_MEAN_RMS_MAX << endm;
126 
127  return StMaker::Init();
128 }
129 
130 //_____________________________________________________________________________
131 void StSvtBadAnodesMaker::setSvtData()
132 {
133  St_DataSet *dataSet;
134 
135  dataSet = GetDataSet("StSvtRawData");
136  assert(dataSet);
137  mSvtData = (StSvtData*)(dataSet->GetObject());
138  assert(mSvtData);
139 }
140 
141 //_____________________________________________________________________________
142 void StSvtBadAnodesMaker::setPedestal()
143 {
144  if (Debug()) gMessMgr->Debug() << "StSvtQAMaker::setPedestal" << endm;
145 
146  St_DataSet *dataSet;
147 
148  dataSet = GetDataSet("StSvtPedestal");
149 
150  if (dataSet) {
151  mSvtPed = (StSvtHybridCollection*)(dataSet->GetObject());
152  assert(mSvtPed);
153  }
154 }
155 
156 //_____________________________________________________________________________
157 void StSvtBadAnodesMaker::setRMSPedestal()
158 {
159  if (Debug()) gMessMgr->Debug() << "StSvtQAMaker::setRMSPedestal" << endm;
160 
161  St_DataSet *dataSet;
162 
163  dataSet = GetDataSet("StSvtRMSPedestal");
164 
165  if (dataSet) {
166  mSvtRMSPed = (StSvtHybridCollection*)(dataSet->GetObject());
167  assert(mSvtRMSPed);
168  }
169 }
170 
171 //_____________________________________________________________________________
172 void StSvtBadAnodesMaker::setSvtBadAnodes()
173 {
174  St_DataSet *dataset;
175 
176  dataset = new TObjectSet("StSvtBadAnodes");
177  AddConst(dataset);
178 
179  if (!mSvtBadAnodes) {
180  mSvtBadAnodes = new StSvtHybridCollection(mSvtData->getConfiguration());
181  dataset->SetObject((TObject*)mSvtBadAnodes);
182  assert(mSvtBadAnodes);
183  }
184 }
185 
186 //_____________________________________________________________________________
187 void StSvtBadAnodesMaker::bookHistograms()
188 {
189  char posTitle[10];
190  char preTitle[25];
191  char* title;
192  int index, index2;
193 
194  mFile = new TFile(mFileName->Data(),"RECREATE");
195 
196  mBadAnodesHist = new TH2F*[mSvtData->getTotalNumberOfHybrids()];
197 
198  mBadAnodesBarrel = new TH1F*[3];
199  mBadAnodesBarrel[0] = new TH1F("barrel1","Fraction of Bad Anodes",8,0.5,8.5);
200  mBadAnodesBarrel[1] = new TH1F("barrel2","Fraction of Bad Anodes",12,0.5,12.5);
201  mBadAnodesBarrel[2] = new TH1F("barrel3","Fraction of Bad Anodes",16,0.5,16.5);
202 
203  mBadAnodesLadder = new TH1F*[36];
204 
205  for (int barrel = 1;barrel <= mSvtData->getNumberOfBarrels();barrel++) {
206  for (int ladder = 1;ladder <= mSvtData->getNumberOfLadders(barrel);ladder++) {
207  for (int wafer = 1;wafer <= mSvtData->getNumberOfWafers(barrel);wafer++) {
208  for (int hybrid = 1;hybrid <= mSvtData->getNumberOfHybrids();hybrid++) {
209 
210  index = mSvtData->getHybridIndex(barrel,ladder,wafer,hybrid);
211  if(index < 0) continue;
212 
213  sprintf(preTitle,"BadAnodes");
214  sprintf(posTitle,"b%dl%dw%dh%d",barrel, ladder, wafer, hybrid);
215  title = strcat(preTitle,posTitle);
216 
217  mBadAnodesHist[index] = new TH2F(title,"Bad Anodes ADC",240,0.5,240.5,128,-0.5,127.5);
218  }
219  }
220 
221  sprintf(preTitle,"ladder");
222  sprintf(posTitle,"b%dl%d",barrel, ladder);
223  title = strcat(preTitle,posTitle);
224 
225  switch (barrel) {
226  case 1:
227  index2 = (ladder-1);
228  mBadAnodesLadder[index2] = new TH1F(title,"Fraction of Bad Anodes",8,0.5,8.5);
229  break;
230  case 2:
231  index2 = mSvtData->getNumberOfLadders(barrel) + (ladder-1);
232  mBadAnodesLadder[index2] = new TH1F(title,"Fraction of Bad Anodes",12,0.5,12.5);
233  break;
234  case 3:
235  index2 = mSvtData->getNumberOfLadders(barrel) + (ladder-1);
236  mBadAnodesLadder[index2] = new TH1F(title,"Fraction of Bad Anodes",14,0.5,14.5);
237  break;
238  }
239  }
240  }
241 
242 }
243 
244 //_____________________________________________________________________________
246 {
247  int adc, anode, nAnodes, nSeq, iseq, time, timeSeq, status, index;
248  int* anodeList;
249  StSequence* Seq;
250  int nullAdc, overloadedAdc, occup, badRMS;
251  float ped,rms,meanPed, meanRms;
252 
253  setSvtData();
254 
255  mEvents++;
256 
257  // scan through the data
258  if (mSvtData) {
259  for (int barrel = 1;barrel <= mSvtData->getNumberOfBarrels();barrel++) {
260  for (int ladder = 1;ladder <= mSvtData->getNumberOfLadders(barrel);ladder++) {
261  for (int wafer = 1;wafer <= mSvtData->getNumberOfWafers(barrel);wafer++) {
262  for (int hybrid = 1;hybrid <= mSvtData->getNumberOfHybrids();hybrid++) {
263 
264  index = mSvtData->getHybridIndex(barrel, ladder, wafer, hybrid);
265  if (index < 0) continue;
266 
267  mHybridData = (StSvtHybridData*)mSvtData->at(index);
268 
269  if (!mHybridData) continue;
270 
271  mHybridBadAnodes = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(index);
272  if (!mHybridBadAnodes)
273  mHybridBadAnodes = new StSvtHybridBadAnodes(barrel, ladder, wafer, hybrid);
274 
275  anodeList = NULL;
276  nAnodes = mHybridData->getAnodeList(anodeList);
277 
278  for (int ianode=0;ianode<nAnodes;ianode++) { // loop over anodes
279 
280  // initialize variables
281  nullAdc = 0;
282  overloadedAdc = 0;
283  occup = 0;
284  badRMS = 0;
285 
286  anode = anodeList[ianode];
287  Seq = NULL;
288  nSeq = 0;
289 
290  status = mHybridData->getSequences(anode,nSeq,Seq);
291 
292  for (iseq=0;iseq<nSeq;iseq++) {
293  for (timeSeq=0; timeSeq<Seq[iseq].length; timeSeq++) {
294 
295  time = Seq[iseq].startTimeBin + timeSeq;
296  adc =(int)Seq[iseq].firstAdc[timeSeq];
297 
298  // check whether adc is NULL
299  if (adc == NULL_ADC)
300  nullAdc++;
301 
302  // check whether adc is overloaded
303  if (adc == OVERLOADED_ADC)
304  overloadedAdc++;
305 
306  // check the occupancy
307  if ((adc > 0) && (adc < 255))
308  occup++;
309  }
310  }
311 
312  if (nullAdc >= NULL_ADC_THRESHOLD)
313  mHybridBadAnodes->addNullAdc(anode,kTRUE,mEvents);
314  else
315  mHybridBadAnodes->addNullAdc(anode,kFALSE,mEvents);
316  if (overloadedAdc >= OVERLOADED_ADC_THRESHOLD)
317  mHybridBadAnodes->addOverloadedAdc(anode,kTRUE,mEvents);
318  else
319  mHybridBadAnodes->addOverloadedAdc(anode,kFALSE,mEvents);
320  if (occup >= OCCUP_THRESHOLD)
321  mHybridBadAnodes->addHighOccup(anode,kTRUE,mEvents);
322  else
323  mHybridBadAnodes->addHighOccup(anode,kFALSE,mEvents);
324  }
325 
326  mSvtBadAnodes->put_at((TObject*)mHybridBadAnodes,index);
327  }
328  }
329  }
330  }
331  }
332 
333  // scan through the data
334  if (mSvtPed) {
335  for (int barrel = 1;barrel <= mSvtPed->getNumberOfBarrels();barrel++) {
336  for (int ladder = 1;ladder <= mSvtPed->getNumberOfLadders(barrel);ladder++) {
337  for (int wafer = 1;wafer <= mSvtPed->getNumberOfWafers(barrel);wafer++) {
338  for (int hybrid = 1;hybrid <= mSvtPed->getNumberOfHybrids();hybrid++) {
339 
340  index = mSvtPed->getHybridIndex(barrel, ladder, wafer, hybrid);
341  if (index < 0) continue;
342 
343  mHybridPed = (StSvtHybridPed*)mSvtPed->at(index);
344 
345  if (!mHybridPed) continue;
346 
347  mHybridBadAnodes = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(index);
348  if (!mHybridBadAnodes)
349  mHybridBadAnodes = new StSvtHybridBadAnodes(barrel, ladder, wafer, hybrid);
350 
351  for (int anode=1;anode<=240;anode++) { // loop over anodes
352  meanPed=0;
353  for (time=4; time<127; time++) {
354 
355  // check pedestal
356  ped = mHybridPed->At(mHybridPed->getPixelIndex(anode, time));
357  meanPed += ped;
358  }
359 
360  meanPed /= 123;
361 
362  if ((meanPed < BAD_MEAN_PED_MIN) || (meanPed > BAD_MEAN_PED_MAX))
363  mHybridBadAnodes->setBadAnode(anode);
364  }
365 
366  mSvtBadAnodes->put_at((TObject*)mHybridBadAnodes,index);
367  }
368  }
369  }
370  }
371  }
372 
373  // scan through the data
374  if (mSvtRMSPed) {
375  for (int barrel = 1;barrel <= mSvtRMSPed->getNumberOfBarrels();barrel++) {
376  for (int ladder = 1;ladder <= mSvtRMSPed->getNumberOfLadders(barrel);ladder++) {
377  for (int wafer = 1;wafer <= mSvtRMSPed->getNumberOfWafers(barrel);wafer++) {
378  for (int hybrid = 1;hybrid <= mSvtRMSPed->getNumberOfHybrids();hybrid++) {
379 
380  index = mSvtRMSPed->getHybridIndex(barrel, ladder, wafer, hybrid);
381  if (index < 0) continue;
382 
383  mHybridRMSPed = (StSvtHybridPed*)mSvtRMSPed->at(index);
384 
385  if (!mHybridRMSPed) continue;
386 
387  mHybridBadAnodes = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(index);
388  if (!mHybridBadAnodes)
389  mHybridBadAnodes = new StSvtHybridBadAnodes(barrel, ladder, wafer, hybrid);
390 
391  for (int anode=1;anode<=240;anode++) { // loop over anodes
392  meanRms=0;
393  for (time=4; time<127; time++) {
394 
395  // check pedestal RMS
396  rms = (mHybridRMSPed->At(mHybridRMSPed->getPixelIndex(anode, time)))/mRmsScaleFactor;
397  meanRms += rms;
398  //if (rms > BAD_RMS)
399  //badRMS++;
400  }
401 
402  meanRms /= 123;
403 
404  //if (badRMS >= RMS_THRESHOLD)
405  //mHybridBadAnodes->setBadAnode(anode);
406  if ((meanRms < BAD_MEAN_RMS_MIN) || (meanRms > BAD_MEAN_RMS_MAX))
407  mHybridBadAnodes->setBadAnode(anode);
408  }
409 
410  mSvtBadAnodes->put_at((TObject*)mHybridBadAnodes,index);
411  }
412  }
413  }
414  }
415  }
416 
417  return kStOK;
418 }
419 
420 //_____________________________________________________________________________
421 void StSvtBadAnodesMaker::writeToFile(const char* fileName)
422 {
423  mFileNameTxt = new TString(mFileName->Data());
424  mFileNameTxt->ReplaceAll("root","txt");
425  ofstream file2(mFileNameTxt->Data());
426  mFileNameTxt->ReplaceAll("txt","Daq.txt");
427  ofstream file(mFileNameTxt->Data());
428 
429  int anode, index, nTotalBadAnodes = 0;
430  int recBoard, mezz, mz_hyb;
431 
432  // scan through the bad anodes data and write them out
433  for (int barrel = 1;barrel <= mSvtBadAnodes->getNumberOfBarrels();barrel++) {
434  for (int ladder = 1;ladder <= mSvtBadAnodes->getNumberOfLadders(barrel);ladder++) {
435  for (int wafer = 1;wafer <= mSvtBadAnodes->getNumberOfWafers(barrel);wafer++) {
436  for (int hybrid = 1;hybrid <= mSvtBadAnodes->getNumberOfHybrids();hybrid++) {
437 
438  blwh2rma(barrel, ladder, wafer, hybrid, recBoard, mezz, mz_hyb);
439 
440  index = mSvtBadAnodes->getHybridIndex(barrel, ladder, wafer, hybrid);
441  if (index < 0) continue;
442 
443  mHybridBadAnodes = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(index);
444 
445  if (!mHybridBadAnodes) continue;
446 
447  for (int ianode=0;ianode<MAX_NUMBER_OF_ANODES;ianode++) {
448  anode = ianode + 1;
449  if (mHybridBadAnodes->isBadAnode(anode)) {
450  file << recBoard << " " << mezz << " " << mz_hyb << " "
451  << anode << " " << "0.0" << endl;
452  file2 << barrel << " " << ladder << " " << wafer << " " << hybrid << " " << anode << endl;
453  nTotalBadAnodes++;
454  }
455  }
456  }
457  }
458  }
459  }
460 
461  file2 << (float)nTotalBadAnodes/((float)MAX_NUMBER_OF_ANODES*432.) << endl;
462 }
463 
464 //_____________________________________________________________________________
466 {
467  if (Debug()) gMessMgr->Debug() << "StSvtBadAnodesMaker::Finish" << endm;
468 
469  cout << FREQ_NULL_ADC << " " << FREQ_OVERLOADED_ADC << " " << FREQ_OCCUP << endl;
470 
471  int adc, anode, nSeq, iseq, time, timeSeq, status, index, index2, ihybrid;
472  int nBadAnodesHyb=0, nBadAnodesLadder=0;
473  StSequence* Seq;
474  float fraction;
475 
476  // scan through the bad anodes data and check frequency of appearance
477  for (int barrel = 1;barrel <= mSvtBadAnodes->getNumberOfBarrels();barrel++) {
478  for (int ladder = 1;ladder <= mSvtBadAnodes->getNumberOfLadders(barrel);ladder++) {
479  for (int wafer = 1;wafer <= mSvtBadAnodes->getNumberOfWafers(barrel);wafer++) {
480  for (int hybrid = 1;hybrid <= mSvtBadAnodes->getNumberOfHybrids();hybrid++) {
481 
482  index = mSvtBadAnodes->getHybridIndex(barrel, ladder, wafer, hybrid);
483  if (index < 0) continue;
484 
485  mHybridBadAnodes = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(index);
486 
487  if (!mHybridBadAnodes) continue;
488 
489  mHybridData = (StSvtHybridData*)mSvtData->at(index);
490 
491  if (!mHybridData) continue;
492 
493  for (int ianode=0;ianode<MAX_NUMBER_OF_ANODES;ianode++) {
494  anode = ianode + 1;
495 
496  if ((mHybridBadAnodes->getNullAdc(anode) >= FREQ_NULL_ADC) ||
497  (mHybridBadAnodes->getOverloadedAdc(anode) >= FREQ_OVERLOADED_ADC) ||
498  (mHybridBadAnodes->getHighOccup(anode) >= FREQ_OCCUP))
499  mHybridBadAnodes->setBadAnode(anode);
500 
501  // fill histograms
502  if (Debug()) {
503 
504  if (mHybridBadAnodes->isBadAnode(anode)) {
505 
506  nBadAnodesHyb++;
507  nBadAnodesLadder++;
508 
509  Seq = NULL;
510  nSeq = 0;
511 
512  status = mHybridData->getSequences(anode,nSeq,Seq);
513 
514  for (iseq=0;iseq<nSeq;iseq++) {
515  for (timeSeq=0; timeSeq<Seq[iseq].length; timeSeq++) {
516 
517  time = Seq[iseq].startTimeBin + timeSeq;
518  adc =(int)Seq[iseq].firstAdc[timeSeq];
519 
520  mBadAnodesHist[index]->Fill(anode,time,adc);
521  }
522  }
523  }
524  }
525  }
526 
527  if (Debug()) {
528  ihybrid = (wafer-1)*2+hybrid;
529  fraction = (float)nBadAnodesHyb/240.;
530 
531  switch (barrel) {
532  case 1:
533  index2 = (ladder-1);
534  mBadAnodesLadder[index2]->Fill(ihybrid,fraction);
535  break;
536  case 2:
537  index2 = mSvtData->getNumberOfLadders(barrel) + (ladder-1);
538  mBadAnodesLadder[index2]->Fill(ihybrid,fraction);
539  break;
540  case 3:
541  index2 = mSvtData->getNumberOfLadders(barrel) + (ladder-1);
542  mBadAnodesLadder[index2]->Fill(ihybrid,fraction);
543  break;
544  }
545 
546  nBadAnodesHyb = 0;
547  }
548  }
549  }
550 
551  if (Debug()) {
552  fraction = (float)nBadAnodesLadder/(240.*(float)mSvtBadAnodes->getNumberOfLadders(barrel));
553  mBadAnodesBarrel[barrel-1]->Fill(ladder,fraction);
554  nBadAnodesLadder = 0;
555  }
556  }
557  }
558 
559  writeToFile();
560 
561  if (Debug()) {
562  mFile->Write();
563  mFile->Close();
564  }
565 
566  return kStOK;
567 }
568 
569 //_____________________________________________________________________________
570 void StSvtBadAnodesMaker::Reset()
571 {
572  if (Debug()) gMessMgr->Debug() << "StSvtBadAnodesMaker::Clear" << endm;
573 
574  int index;
575  for (int barrel = 1;barrel <= mSvtBadAnodes->getNumberOfBarrels();barrel++) {
576  for (int ladder = 1;ladder <= mSvtBadAnodes->getNumberOfLadders(barrel);ladder++) {
577  for (int wafer = 1;wafer <= mSvtBadAnodes->getNumberOfWafers(barrel);wafer++) {
578  for (int hybrid = 1;hybrid <= mSvtBadAnodes->getNumberOfHybrids();hybrid++) {
579 
580  index = mSvtBadAnodes->getHybridIndex(barrel, ladder, wafer, hybrid);
581  if (index < 0) continue;
582 
583  mHybridBadAnodes = (StSvtHybridBadAnodes*)mSvtBadAnodes->at(index);
584 
585  if (mHybridBadAnodes) delete mHybridBadAnodes;
586  mSvtBadAnodes->put_at(NULL,index);
587  }
588  }
589  }
590  }
591 }
592 
593 //_____________________________________________________________________________
594 void StSvtBadAnodesMaker::blwh2rma(int barrel, int ladder, int wafer, int hybrid,
595  int& recBoard, int& mezz, int& mz_hyb)
596 {
597  switch (barrel) {
598  case 1:
599  recBoard = key_barrel1[ladder-1][wafer-1][hybrid-1][0];
600  mezz = key_barrel1[ladder-1][wafer-1][hybrid-1][1];
601  mz_hyb = key_barrel1[ladder-1][wafer-1][hybrid-1][2];
602  break;
603  case 2:
604  recBoard = key_barrel2[ladder-1][wafer-1][hybrid-1][0];
605  mezz = key_barrel2[ladder-1][wafer-1][hybrid-1][1];
606  mz_hyb = key_barrel2[ladder-1][wafer-1][hybrid-1][2];
607  break;
608  case 3:
609  recBoard = key_barrel3[ladder-1][wafer-1][hybrid-1][0];
610  mezz = key_barrel3[ladder-1][wafer-1][hybrid-1][1];
611  mz_hyb = key_barrel3[ladder-1][wafer-1][hybrid-1][2];
612  break;
613  }
614 }
virtual void SetObject(TObject *obj)
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:480
virtual TObject * GetObject() const
The depricated method (left here for the sake of the backward compatibility)
Definition: TDataSet.cxx:428
Definition: Stypes.h:40