]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/Alieve/MUONData.cxx
Added interface to load-time threshold/(auto)pedestal settings.
[u/mrichter/AliRoot.git] / EVE / Alieve / MUONData.cxx
1 //
2 // Sources:
3 //
4 // GetTrackerMapping = AliMUONDigitMaker::GetMapping
5 // GetTriggerMapping = AliMUONDigitMaker::TriggerDigits
6 // GetTriggerChamber = AliMUONDigitMaker::GetTriggerChamber
7 // LoadRawTracker    = MUONRawStreamTracker.C
8 // LoadRawTrigger    = MUONRawStreamTrigger.C
9 //
10
11 #include "MUONData.h"
12
13 #include <Alieve/MUONChamberData.h>
14 #include <Alieve/EventAlieve.h>
15
16 #include <AliRawReader.h>
17 #include <AliRawReaderFile.h>
18 #include <AliRawReaderDate.h>
19 #include <AliRawReaderRoot.h>
20
21 #include <AliTracker.h>
22 #include <AliMagFMaps.h>
23
24 #include <AliMUONTrack.h>
25 #include <AliMUONTrackParam.h>
26 #include <AliMUONDigit.h>
27 #include <AliMUONRawStreamTracker.h>
28 #include <AliMUONRawStreamTrigger.h>
29 #include <AliMUONDDLTracker.h>
30 #include <AliMUONBlockHeader.h>
31 #include <AliMUONDspHeader.h>
32 #include <AliMUONBusStruct.h>
33 #include <AliMUONDDLTrigger.h>
34 #include <AliMUONDarcHeader.h>
35 #include <AliMUONRegHeader.h>
36 #include <AliMUONLocalStruct.h>
37 #include <AliMUONTriggerCrateStore.h>
38 #include <AliMUONTriggerCrate.h>
39 #include <AliMUONLocalTriggerBoard.h>
40 #include <AliMUONTriggerCircuit.h>
41 #include <mapping/AliMpBusPatch.h>
42 #include <mapping/AliMpVSegmentation.h>
43 #include <mapping/AliMpSegmentation.h>
44 #include <mapping/AliMpPad.h>
45 #include <mapping/AliMpDEManager.h>
46
47 #include "TTree.h"
48 #include "TString.h"
49 #include "TMatrixD.h"
50 #include "TClonesArray.h"
51
52 using namespace Reve;
53 using namespace Alieve;
54
55 //______________________________________________________________________
56 // MUONData
57 //
58
59 ClassImp(MUONData)
60
61 AliRawReader*            MUONData::fgRawReader        = 0;
62 AliMUONRawStreamTracker* MUONData::fgRawStreamTracker = 0;
63 AliMUONRawStreamTrigger* MUONData::fgRawStreamTrigger = 0;
64 AliMpBusPatch*           MUONData::fgBusPatchManager  = 0;
65
66 //______________________________________________________________________
67 MUONData::MUONData() :
68   fChambers(14),
69   fNTracks(0),
70   fTrackPoints(0),
71   fNPoints(0)
72 {
73   //
74   // Constructor
75   //
76   
77   CreateAllChambers();
78
79 }
80
81 //______________________________________________________________________
82 MUONData::~MUONData()
83 {
84   //
85   // Destructor
86   //
87
88   DeleteAllChambers();
89
90   delete [] fTrackPoints;
91
92   fTrackPoints = 0;
93
94 }
95
96 //______________________________________________________________________
97 MUONData::MUONData(const MUONData &mdata) :
98   TObject(mdata),
99   Reve::ReferenceCount()
100 {
101   //
102   // Copy constructor
103   //
104
105 }
106
107 //______________________________________________________________________
108 MUONData& MUONData::operator=(const MUONData &mdata)
109 {
110   //
111   // Assignment operator
112   //
113
114   if (this != &mdata) {
115
116   }
117
118   return *this;
119
120 }
121
122 //______________________________________________________________________
123 void MUONData::CreateChamber(Int_t chamber)
124 {
125   // 
126   // create data for the chamber with id=chamber (0 to 13)
127   //
128
129   if (fChambers[chamber] == 0)
130     fChambers[chamber] = new MUONChamberData(chamber);
131
132 }
133
134 //______________________________________________________________________
135 void MUONData::CreateAllChambers()
136 {
137   //
138   // create all 14 chambers data
139   //
140
141   for (Int_t c = 0; c < 14; ++c)
142     CreateChamber(c);
143
144 }
145
146 //______________________________________________________________________
147 void MUONData::DropAllChambers()
148 {
149   // 
150   // release data from all chambers 
151   //
152
153   for (Int_t c = 0; c < 14; ++c) {
154
155     if (fChambers[c] != 0)
156       fChambers[c]->DropData();
157
158   }
159
160 }
161
162 //______________________________________________________________________
163 void MUONData::DeleteAllChambers()
164 {
165   //
166   // delete all chambers data
167   //
168
169   for (Int_t c = 0; c < 14; ++c) {
170
171     delete fChambers[c];
172     fChambers[c] = 0;
173
174   }
175
176 }
177
178 //______________________________________________________________________
179 void MUONData::LoadTracks(TTree* tree)
180 {
181   //
182   // load tracks from the TreeT and creates the track points array
183   // the structure of fTrackPoints:
184   // 0,  0,  0, - new track
185   // px, py, pz - from track param at vertex
186   // x, y, z    - track param at vertex
187   // x, y, z    - track param at hits
188   // ..........
189   //
190   // 0,  0,  0  - new track
191   // ..........
192   //
193
194   Int_t maxTrackHits = 1+2*10+4;
195
196   TClonesArray *tracks = 0;
197   tree->SetBranchAddress("MUONTrack",&tracks);
198   tree->GetEntry(0);
199
200   Int_t ntracks = tracks->GetEntriesFast();
201   printf("Found %d tracks. \n",ntracks);
202
203   fNTracks = ntracks;
204
205   Int_t maxTrackPoints = (3+3+3*maxTrackHits)*ntracks;
206   fTrackPoints = new Float_t [maxTrackPoints];
207
208   TMatrixD smatrix(2,2);
209   TMatrixD sums(2,1);
210   TMatrixD res(2,1);
211
212   Float_t xRec, xRec0;
213   Float_t yRec, yRec0;
214   Float_t zRec, zRec0;
215   Float_t px0, py0, pz0;
216   
217   Float_t zg[4] = { -1603.5, -1620.5, -1703.5, -1720.5 };
218
219   AliMUONTrack *mt;  
220   Int_t count = 0;
221   for (Int_t n = 0; n < ntracks; n++) {
222
223     if (count >= maxTrackPoints) continue;
224     fTrackPoints[3*count  ] = 0.0;
225     fTrackPoints[3*count+1] = 0.0;
226     fTrackPoints[3*count+2] = 0.0;
227     count++;
228
229     mt = (AliMUONTrack*) tracks->At(n);
230
231     printf("Match trigger %d \n",mt->GetMatchTrigger());
232
233     AliMUONTrackParam *trackParam = mt->GetTrackParamAtVertex(); 
234     xRec0  = trackParam->GetNonBendingCoor();
235     yRec0  = trackParam->GetBendingCoor();
236     zRec0  = trackParam->GetZ();
237
238     px0 = trackParam->Px();
239     py0 = trackParam->Py();
240     pz0 = trackParam->Pz();
241
242     if (count >= maxTrackPoints) continue;
243     fTrackPoints[3*count  ] = px0;
244     fTrackPoints[3*count+1] = py0;
245     fTrackPoints[3*count+2] = pz0;
246     count++;
247
248     if (count >= maxTrackPoints) continue;
249     fTrackPoints[3*count  ] = xRec0;
250     fTrackPoints[3*count+1] = yRec0;
251     fTrackPoints[3*count+2] = zRec0;
252     count++;
253     
254     Float_t xr[20], yr[20], zr[20];
255     for (Int_t i = 0; i < 10; i++) xr[i]=yr[i]=zr[i]=0.0;
256
257     Int_t nTrackHits = mt->GetNTrackHits();
258     printf("Nhits = %d \n",nTrackHits);
259     TClonesArray* trackParamAtHit;
260     for (Int_t iHit = 0; iHit < nTrackHits; iHit++){
261       trackParamAtHit = mt->GetTrackParamAtHit();
262       trackParam = (AliMUONTrackParam*) trackParamAtHit->At(iHit); 
263       xRec  = trackParam->GetNonBendingCoor();
264       yRec  = trackParam->GetBendingCoor();
265       zRec  = trackParam->GetZ();
266
267       //printf("Hit %d x %f y %f z %f \n",iHit,xRec,yRec,zRec);
268
269       xr[iHit] = xRec;
270       yr[iHit] = yRec;
271       zr[iHit] = zRec;
272
273       if (count >= maxTrackPoints) continue;
274       fTrackPoints[3*count  ] = xRec;
275       fTrackPoints[3*count+1] = yRec;
276       fTrackPoints[3*count+2] = zRec;
277       count++;
278     
279     }
280
281     Float_t xrc[20], yrc[20], zrc[20];
282     Int_t nrc = 0;
283     if (mt->GetMatchTrigger() && 1) {
284
285       for (Int_t i = 0; i < nTrackHits; i++) {
286         if (TMath::Abs(zr[i]) > 1000.0) {
287           //printf("Hit %d x %f y %f z %f \n",iHit,xr[i],yr[i],zr[i]);
288           xrc[nrc] = xr[i];
289           yrc[nrc] = yr[i];
290           zrc[nrc] = zr[i];
291           nrc++;
292         }
293       }
294
295       if (nrc < 2) continue;
296
297       Double_t xv, yv;
298       Float_t ax, bx, ay, by;
299       
300       // fit x-z
301       smatrix.Zero();
302       sums.Zero();
303       for (Int_t i = 0; i < nrc; i++) {
304         xv = (Double_t)zrc[i];
305         yv = (Double_t)xrc[i];
306         //printf("x-z: xv %f yv %f \n",xv,yv);
307         smatrix(0,0) += 1.0;
308         smatrix(1,1) += xv*xv;
309         smatrix(0,1) += xv;
310         smatrix(1,0) += xv;
311         sums(0,0)    += yv;
312         sums(1,0)    += xv*yv;
313       }
314       res = smatrix.Invert() * sums;
315       ax = res(0,0);
316       bx = res(1,0);
317
318       // fit y-z
319       smatrix.Zero();
320       sums.Zero();
321       for (Int_t i = 0; i < nrc; i++) {
322         xv = (Double_t)zrc[i];
323         yv = (Double_t)yrc[i];
324         //printf("y-z: xv %f yv %f \n",xv,yv);
325         smatrix(0,0) += 1.0;
326         smatrix(1,1) += xv*xv;
327         smatrix(0,1) += xv;
328         smatrix(1,0) += xv;
329         sums(0,0)    += yv;
330         sums(1,0)    += xv*yv;
331       }
332       res = smatrix.Invert() * sums;
333       ay = res(0,0);
334       by = res(1,0);
335
336       Float_t xtc, ytc, ztc;
337       for (Int_t ii = 0; ii < 4; ii++) {
338
339         ztc = zg[ii];
340         ytc = ay+by*zg[ii];
341         xtc = ax+bx*zg[ii];
342
343         //printf("tc: x %f y %f z %f \n",xtc,ytc,ztc);
344
345         if (count >= maxTrackPoints) continue;
346         fTrackPoints[3*count  ] = xtc;
347         fTrackPoints[3*count+1] = ytc;
348         fTrackPoints[3*count+2] = ztc;
349         count++;
350
351       }
352
353     }  // end match trigger
354
355   }
356
357   fNPoints = 3*count;
358
359   printf("MUONData found %d track points. \n",fNPoints);
360
361 }
362
363 //______________________________________________________________________
364 void MUONData::LoadDigits(TTree* tree)
365 {
366   // 
367   // load digits from the TreeD
368   //
369
370   Char_t branchname[30];
371   TClonesArray *digits = 0;
372   Int_t ndigits;
373   AliMUONDigit  *mdig;
374   Int_t cathode, detElemId, ix, iy, charge;
375
376   for (Int_t c = 0; c < 14; ++c) {
377
378     if (fChambers[c] == 0) continue;
379     sprintf(branchname,"MUONDigits%d",c+1);
380     tree->SetBranchAddress(branchname,&digits);
381     tree->GetEntry(0);
382
383     ndigits = digits->GetEntriesFast(); 
384
385     for (Int_t id = 0; id < ndigits; id++) {
386       mdig  = (AliMUONDigit*)digits->UncheckedAt(id);
387
388       cathode   = mdig->Cathode();
389       ix        = mdig->PadX();
390       iy        = mdig->PadY();
391       detElemId = mdig->DetElemId();      
392       charge    = (Int_t)mdig->Signal();
393
394       if (c > 9) {
395         //printf("cha %d deid %d cath %1d ix %d iy %d q %d \n",c,detElemId,cathode,ix,iy,charge);  
396       }
397
398       fChambers[c]->RegisterDigit(detElemId,cathode,ix,iy,charge);
399
400     } // end digits loop
401
402   }
403
404 }
405
406 //______________________________________________________________________
407 void MUONData::LoadRaw(TString fileName)
408 {
409   //
410   // load raw data from fileName; tracker and trigger data
411   //
412
413   if (fgRawReader == 0) {
414     // check extention to choose the rawdata file format
415     if (fileName.EndsWith("/")) {
416       fgRawReader = new AliRawReaderFile(fileName); // DDL files
417     } else if (fileName.EndsWith(".root")) {
418       fgRawReader = new AliRawReaderRoot(fileName); // ROOT file
419     } else if (!fileName.IsNull()) {
420       fgRawReader = new AliRawReaderDate(fileName); // DATE file
421     }
422     fgRawStreamTracker = new AliMUONRawStreamTracker(fgRawReader);
423     fgRawStreamTrigger = new AliMUONRawStreamTrigger(fgRawReader);
424     fgBusPatchManager = new AliMpBusPatch();
425     fgBusPatchManager->ReadBusPatchFile();
426   }
427   
428   LoadRawTracker();
429   LoadRawTrigger();
430
431 }
432
433 //______________________________________________________________________
434 void MUONData::LoadRawTracker()
435 {
436   //
437   // load raw data for the tracking chambers
438   //
439
440   fgRawReader->RewindEvents();
441
442   AliMUONDigit* digit = new AliMUONDigit();
443
444   Int_t maxEvent = 1000;
445   Int_t minDDL = 0, maxDDL = 19;
446   Int_t cathode, detElemId, ix, iy, iChamber;
447
448   AliMUONDDLTracker*       ddlTracker = 0x0;
449   AliMUONBlockHeader*      blkHeader  = 0x0;
450   AliMUONDspHeader*        dspHeader  = 0x0;
451   AliMUONBusStruct*        busStruct  = 0x0;
452
453   Int_t iEvent = 0;
454   Int_t dataSize, buspatchId;
455   
456   Event* aevent = Alieve::gEvent;
457
458   while (fgRawReader->NextEvent()) {
459     
460     if (iEvent != aevent->GetEventId()) {
461       iEvent++;
462       continue;
463     }
464     
465     if (iEvent == maxEvent)
466       break;
467     
468     // read DDL while < 20 DDL
469     while(fgRawStreamTracker->NextDDL()) {
470       
471       if (fgRawStreamTracker->GetDDL() < minDDL || 
472           fgRawStreamTracker->GetDDL() > maxDDL)
473         continue;
474       
475       //printf("\niDDL %d\n", fgRawStreamTracker->GetDDL());
476       
477       ddlTracker =  fgRawStreamTracker->GetDDLTracker();
478       
479       // loop over block structure
480       Int_t nBlock = ddlTracker->GetBlkHeaderEntries();
481       for(Int_t iBlock = 0; iBlock < nBlock ;iBlock++){
482         
483         blkHeader = ddlTracker->GetBlkHeaderEntry(iBlock);
484         //printf("Block Total length %d\n",blkHeader->GetTotalLength());
485         
486         // loop over DSP structure
487         Int_t nDsp = blkHeader->GetDspHeaderEntries();
488         for(Int_t iDsp = 0; iDsp < nDsp ;iDsp++){   //DSP loop
489           
490           dspHeader =  blkHeader->GetDspHeaderEntry(iDsp);
491           //   printf("Dsp length %d even word %d\n",dspHeader->GetTotalLength(), dspHeader->GetEventWord());
492           
493           // loop over BusPatch structure
494           Int_t nBusPatch = dspHeader->GetBusPatchEntries();
495           for(Int_t iBusPatch = 0; iBusPatch < nBusPatch; iBusPatch++) {  
496             
497             busStruct = dspHeader->GetBusPatchEntry(iBusPatch);
498             
499             //printf("busPatchId %d", busStruct->GetBusPatchId());
500             //printf(" BlockId %d", busStruct->GetBlockId());
501             //printf(" DspId %d\n", busStruct->GetDspId());
502             
503             // loop over data
504             dataSize = busStruct->GetLength();
505             buspatchId = busStruct->GetBusPatchId();
506             for (Int_t iData = 0; iData < dataSize; iData++) {
507               
508               Int_t  manuId    = busStruct->GetManuId(iData);
509               Int_t  channelId = busStruct->GetChannelId(iData);
510               Int_t  charge    = busStruct->GetCharge(iData);
511               //printf("manuId: %d, channelId: %d charge: %d\n", manuId, channelId, charge);
512               // set digit charge
513               digit->SetSignal(charge);
514               digit->SetPhysicsSignal(charge);
515               digit->SetADC(charge);
516               // Get Back the hits at pads
517               Int_t error;
518               error = GetTrackerMapping(buspatchId,manuId,channelId,digit); 
519               if (error) {
520                 printf("Mapping Error\n");
521                 continue;
522               }
523
524               cathode   = digit->Cathode();
525               ix        = digit->PadX();
526               iy        = digit->PadY();
527               detElemId = digit->DetElemId();      
528               charge    = (Int_t)digit->Signal();
529               iChamber  = detElemId/100 - 1;
530
531               fChambers[iChamber]->RegisterDigit(detElemId,cathode,ix,iy,charge);
532               
533             } // iData
534           } // iBusPatch
535         } // iDsp
536       } // iBlock
537     } // NextDDL
538
539     break;
540
541   }  // end event loop
542   
543   delete digit;
544
545 }
546
547 //______________________________________________________________________
548 void MUONData::LoadRawTrigger()
549 {
550   // 
551   // load raw data for the trigger chambers
552   //
553
554   fgRawReader->RewindEvents();
555
556   Int_t maxEvent = 1000;
557   Int_t minDDL = 0, maxDDL = 1;
558   Int_t detElemId, iChamber, cathode, charge, ix, iy;
559
560   AliMUONDDLTrigger*       ddlTrigger  = 0x0;
561   AliMUONDarcHeader*       darcHeader  = 0x0;
562   AliMUONRegHeader*        regHeader   = 0x0;
563   AliMUONLocalStruct*      localStruct = 0x0;
564   
565   // crate manager
566   AliMUONTriggerCrateStore* crateManager = new AliMUONTriggerCrateStore();   
567   crateManager->ReadFromFile();
568
569   // Loop over events  
570   Int_t iEvent = 0;
571   TList digitList;
572   
573   Event* aevent = Alieve::gEvent;
574   
575   while (fgRawReader->NextEvent()) {
576     
577     if (iEvent != aevent->GetEventId()) {
578       iEvent++;
579       continue;
580     }
581     
582     if (iEvent == maxEvent)
583       break;
584     
585     // read DDL while < 2 DDL
586     while(fgRawStreamTrigger->NextDDL()) {
587       
588       if (fgRawStreamTrigger->GetDDL() < minDDL || 
589           fgRawStreamTrigger->GetDDL() > maxDDL)
590         continue;
591       
592       //printf("\niDDL %d\n", fgRawStreamTrigger->GetDDL());
593       
594       ddlTrigger = fgRawStreamTrigger->GetDDLTrigger();
595       darcHeader = ddlTrigger->GetDarcHeader();
596       
597       //printf("Global output %x\n", (Int_t)darcHeader->GetGlobalOutput());
598       
599       // loop over regional structures
600       Int_t nReg = darcHeader->GetRegHeaderEntries();
601       for(Int_t iReg = 0; iReg < nReg ;iReg++){   //REG loop
602         
603         //printf("RegionalId %d\n", iReg);
604         
605         regHeader =  darcHeader->GetRegHeaderEntry(iReg);
606         //  printf("Reg length %d\n",regHeader->GetHeaderLength());
607         
608         // crate info
609         AliMUONTriggerCrate* crate = crateManager->Crate(fgRawStreamTrigger->GetDDL(), iReg);
610         TObjArray *boards = crate->Boards();
611         
612         // loop over local structures
613         Int_t nLocal = regHeader->GetLocalEntries();
614         for(Int_t iLocal = 0; iLocal < nLocal; iLocal++) {  
615           
616           localStruct = regHeader->GetLocalEntry(iLocal);
617           
618           // check if trigger 
619           if (localStruct->GetTriggerY() == 0) { // no empty data
620             
621             // local trigger circuit number
622             AliMUONLocalTriggerBoard* localBoard = (AliMUONLocalTriggerBoard*)boards->At(iLocal+1);
623             
624             //printf("LocalId %d\n", localStruct->GetId());
625             /*
626             Int_t iLocCard  = localBoard->GetNumber();
627             Int_t loStripX  = (Int_t)localStruct->GetXPos();
628             Int_t loStripY  = (Int_t)localStruct->GetYPos();
629             Int_t loDev     = (Int_t)localStruct->GetXDev();
630             */
631             //printf("iLocCard: %d, XPos: %d, YPos: %d Dev: %d\n", iLocCard, loStripX, loStripY, loDev);
632
633             digitList.Clear();
634             if ( GetTriggerMapping(localBoard, localStruct, digitList) ) {
635               for (Int_t iEntry = 0; iEntry < digitList.GetEntries(); iEntry++) {
636
637                 AliMUONDigit* digit = (AliMUONDigit*)digitList.At(iEntry);
638                 cathode   = digit->Cathode();
639                 ix        = digit->PadX();
640                 iy        = digit->PadY();
641                 detElemId = digit->DetElemId();      
642                 charge    = (Int_t)digit->Signal();
643                 iChamber  = detElemId/100 - 1;
644                 
645                 //printf("cha %d deid %d cath %1d ix %d iy %d q %d \n",iChamber,detElemId,cathode,ix,iy,charge);  
646
647                 fChambers[iChamber]->RegisterDigit(detElemId,cathode,ix,iy,charge);
648
649               }
650
651             }
652
653           }
654         } // iLocal
655       } // iReg
656     } // NextDDL
657
658     break;
659
660   }  // end event loop
661
662   delete crateManager;
663
664 }
665
666 //______________________________________________________________________
667 Int_t MUONData::GetTrackerMapping(Int_t buspatchId, UShort_t manuId, UChar_t channelId, AliMUONDigit* digit)
668 {
669   //
670   // decode digits mapping for the tracking chambers
671   //
672   
673   // getting DE from buspatch
674   Int_t detElemId = fgBusPatchManager->GetDEfromBus(buspatchId);
675   //AliDebug(3,Form("detElemId: %d busPatchId %d\n", detElemId, buspatchId));
676
677   const AliMpVSegmentation* seg = AliMpSegmentation::Instance()->GetMpSegmentationByElectronics(detElemId, manuId);  
678   AliMpPad pad = seg->PadByLocation(AliMpIntPair(manuId,channelId),kTRUE);
679
680   if (!pad.IsValid())
681   {
682     printf("No pad for detElemId: %d, busPatchId %d, manuId: %d, channelId: %d\n",detElemId, buspatchId, manuId, channelId);
683     
684     return 1;
685   } // return error
686
687   // Getting padX, padY and cathode number.
688   Int_t padX = pad.GetIndices().GetFirst();
689   Int_t padY = pad.GetIndices().GetSecond();
690   Int_t iCath = AliMpDEManager::GetCathod(detElemId,seg->PlaneType());
691
692   // storing into digits
693   digit->SetPadX(padX);
694   digit->SetPadY(padY);
695   digit->SetCathode(iCath);
696   digit->SetDetElemId(detElemId);
697   digit->SetElectronics(manuId,channelId);
698   
699   //printf("detElemId: %d, busPatchId %d, manuId: %d, channelId: %d, padx: %d pady %d\n",detElemId, buspatchId, manuId, channelId, padX, padY);
700   
701   return 0;
702
703 }
704
705 //______________________________________________________________________
706 Int_t MUONData::GetTriggerMapping(AliMUONLocalTriggerBoard* localBoard, 
707                                   AliMUONLocalStruct* localStruct,
708                                   TList& digitList)
709 {
710   //
711   // decode digits mapping for the trigger chambers
712   //
713
714   Int_t detElemId;
715   Int_t nBoard;
716   Int_t iCath = -1;
717   Int_t iChamber = 0;
718   Int_t xyPattern = 0;
719
720   // loop over x1-4 and y1-4
721   for (Int_t icase = 0; icase < 8; icase++) {
722
723     // get chamber, cathode and associated trigger response pattern
724     GetTriggerChamber(localStruct, xyPattern, iChamber, iCath, icase);
725   
726     if (!xyPattern) continue;
727
728     // get detElemId
729     AliMUONTriggerCircuit triggerCircuit;
730     detElemId = triggerCircuit.DetElemId(iChamber, localBoard->GetName());
731     nBoard    = localBoard->GetNumber();
732
733     const AliMpVSegmentation* seg = AliMpSegmentation::Instance()->GetMpSegmentation(detElemId, iCath);  
734
735     // loop over the 16 bits of pattern
736     for (Int_t ibitxy = 0; ibitxy < 16; ibitxy++) {
737     
738       if ((xyPattern >> ibitxy) & 0x1) {
739
740         // not quite sure about this
741         Int_t offset = 0;
742         if (iCath && localBoard->GetSwitch(6)) offset = -8;
743
744         AliMpPad pad = seg->PadByLocation(AliMpIntPair(nBoard,ibitxy+offset),kTRUE);
745
746         AliMUONDigit* digit = new  AliMUONDigit();
747         if (!pad.IsValid()) {
748           AliWarning(Form("No pad for detElemId: %d, nboard %d, ibitxy: %d\n",
749                           detElemId, nBoard, ibitxy));
750           continue;
751         } // 
752
753         Int_t padX = pad.GetIndices().GetFirst();
754         Int_t padY = pad.GetIndices().GetSecond();
755
756         // file digit
757         digit->SetSignal(1);
758         digit->SetPadX(padX);
759         digit->SetPadY(padY);
760         digit->SetCathode(iCath);
761         digit->SetDetElemId(detElemId);
762         digit->SetElectronics(nBoard, ibitxy);
763         digitList.Add(digit);
764         
765       }// xyPattern
766     }// ibitxy
767   }// case
768
769   return 1;
770
771 }
772
773 //____________________________________________________________________
774 void MUONData::GetTriggerChamber(AliMUONLocalStruct* localStruct, Int_t& xyPattern, Int_t& iChamber, Int_t& iCath, Int_t icase)
775 {
776   //
777   // extract digits pattern
778   //  
779
780   // get chamber & cathode number, (chamber starts at 0 !)
781   switch(icase) {
782   case 0: 
783     xyPattern =  localStruct->GetX1();
784     iCath = 0;
785     iChamber = 10;
786     break;
787   case 1: 
788     xyPattern =  localStruct->GetX2();
789     iCath = 0;
790     iChamber = 11;
791     break;
792   case 2: 
793     xyPattern =  localStruct->GetX3();
794     iCath = 0;
795     iChamber = 12;
796     break;
797   case 3: 
798     xyPattern =  localStruct->GetX4();
799     iCath = 0;
800     iChamber = 13;
801     break;
802   case 4: 
803     xyPattern =  localStruct->GetY1();
804     iCath = 1;
805     iChamber = 10;
806     break;
807   case 5: 
808     xyPattern =  localStruct->GetY2();
809     iCath = 1;
810     iChamber = 11;
811     break;
812   case 6: 
813     xyPattern =  localStruct->GetY3();
814     iCath = 1;
815     iChamber = 12;
816     break;
817   case 7: 
818     xyPattern =  localStruct->GetY4();
819     iCath = 1;
820     iChamber = 13;
821     break;
822   }
823
824 }
825
826 //______________________________________________________________________
827 MUONChamberData* MUONData::GetChamberData(Int_t chamber)
828 {
829   //
830   // return chamber data
831   //
832
833   if (chamber < 0 || chamber > 13) return 0;
834
835   //if (fChambers[chamber] == 0) CreateChamber(chamber);
836
837   return fChambers[chamber];
838
839 }