]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EVE/EveDet/AliEveTOFDigitsInfo.cxx
coverity
[u/mrichter/AliRoot.git] / EVE / EveDet / AliEveTOFDigitsInfo.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2008, ALICE Experiment at CERN, all rights reserved. *
3  * See http://aliceinfo.cern.ch/Offline/AliRoot/License.html for          *
4  * full copyright notice.                                                 *
5  **************************************************************************/
6
7 // $Id$
8 // Main authors: Matevz Tadel & Alja Mrak-Tadel: 2006, 2007
9
10 //
11 // Class to map TOF digit/raw data information
12 //
13 // Author: A. De Caro (email: decaro@sa.infn.it)
14 //
15
16 #include <TClonesArray.h>
17 #include <TTree.h>
18
19 //#include <TEveTreeTools.h>
20
21 #include <AliDAQ.h>
22 #include <AliLog.h>
23 #include <AliRawReader.h>
24
25 #include <AliTOFCableLengthMap.h>
26 #include <AliTOFdigit.h>
27 #include <AliTOFGeometry.h>
28 #include <AliTOFrawData.h>
29 #include <AliTOFRawStream.h>
30 #include <AliTOFDigitMap.h>
31
32 #include "AliEveTOFDigitsInfo.h"
33
34 //_________________________________________________________
35
36 ClassImp(AliEveTOFDigitsInfo)
37
38   AliEveTOFDigitsInfo::AliEveTOFDigitsInfo(): 
39     TObject(),
40     TEveRefCnt(),
41     fTree (0),
42     fNewTree (0),
43     fGeom (new AliTOFGeometry()),
44     fTOFdigitMap(new AliTOFDigitMap())
45 {}
46 /* ******************************************************* */
47
48 AliEveTOFDigitsInfo:: ~AliEveTOFDigitsInfo() 
49 {
50   //dtr
51
52   delete fGeom;
53   delete fTree;
54   delete fNewTree;
55   delete fTOFdigitMap;
56
57 }
58 /* ******************************************************* */
59
60 void AliEveTOFDigitsInfo::SetTree(TTree * const tree)
61 {
62   //
63   // Set fTree global variable
64   //
65
66   static const TEveException kEH("AliEveTOFDigitsInfo::SetTree ");
67   
68   if(fGeom == 0) {
69     fGeom = new AliTOFGeometry();
70   }
71
72   fTree = tree;
73   /*
74   DecRefCount();
75   IncRefCount();
76   */
77 }
78 /* ******************************************************* */
79 void AliEveTOFDigitsInfo::ReadRaw(AliRawReader* rawReader, Int_t newDecoder)
80 {
81   //
82   // Read raw-data. AliTOFdigit is used to
83   // store raw-adata for all sub-detectors.
84   //
85
86   //AliTOFCableLengthMap *cableLength = new AliTOFCableLengthMap();
87
88   //ofstream ftxt;
89   //Char_t fileName[100];
90   //sprintf(fileName,"TOFrawDataReadingFromEVE%d.txt",nEvent);
91
92   //ftxt.open(fileName,ios::app);
93   //ftxt << endl;
94   //ftxt << "  " << nEvent << endl;
95
96   //if (nEvent<0) printf("%3i\n", nEvent); // only to use nEvent variable
97
98   const Int_t kDDL = AliDAQ::NumberOfDdls("TOF");
99
100   TClonesArray *tofDigits = new TClonesArray("AliTOFdigit",10000);
101   fTree = new TTree();
102   fTree->Branch("TOF", &tofDigits, 32000);
103   fTree->GetEntry(0);
104
105   TClonesArray * clonesRawData = 0x0;
106
107   Int_t detectorIndex[5];
108   Int_t digit[4];
109
110   AliTOFRawStream stream(rawReader);
111
112   for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++) {
113
114     rawReader->Reset();
115     if (newDecoder==0) stream.LoadRawData(indexDDL);
116     else if (newDecoder==1) stream.LoadRawDataBuffers(indexDDL);
117     else if (newDecoder==2) stream.LoadRawDataBuffersV2(indexDDL);
118
119     clonesRawData = (TClonesArray*)stream.GetRawData();
120
121     if (clonesRawData->GetEntriesFast()) AliDebug(2, Form(" Number of TOF digits in the sector number %2i: %5i", indexDDL, clonesRawData->GetEntriesFast()));
122
123     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
124
125       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
126
127       if (tofRawDatum->GetTOF()==-1) continue;
128
129       //Int_t cLenInt = Int_t(cableLength->GetCableTimeShift(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),tofRawDatum->GetTDC())*1000./AliTOFGeometry::TdcBinWidth());
130       digit[0] = tofRawDatum->GetTOF();// - cLenInt;
131       digit[1] = tofRawDatum->GetTOT();
132       digit[2] = tofRawDatum->GetTOT();
133       digit[3] = -1;
134
135       /*
136       if (indexDDL<10) ftxt << "  " << indexDDL;
137       else             ftxt << " " << indexDDL;
138       if (tofRawDatum->GetTRM()<10) ftxt << "  " << tofRawDatum->GetTRM();
139       else                          ftxt << " " << tofRawDatum->GetTRM();
140       ftxt << "  " << tofRawDatum->GetTRMchain();
141       if (tofRawDatum->GetTDC()<10) ftxt << "  " << tofRawDatum->GetTDC();
142       else                          ftxt << " " << tofRawDatum->GetTDC();
143       ftxt << "  " << tofRawDatum->GetTDCchannel();
144       */
145
146       stream.EquipmentId2VolumeId(indexDDL, tofRawDatum->GetTRM(), tofRawDatum->GetTRMchain(),
147                                   tofRawDatum->GetTDC(), tofRawDatum->GetTDCchannel(), detectorIndex);
148
149       /* check valid index */
150       if (detectorIndex[0]==-1||detectorIndex[1]==-1||detectorIndex[2]==-1||detectorIndex[3]==-1||detectorIndex[4]==-1) continue;
151
152       // Do not reconstruct anything in the holes
153       if (detectorIndex[0]==13 || detectorIndex[0]==14 || detectorIndex[0]==15 ) { // sectors with holes
154         if (detectorIndex[1]==2) { // plate with holes
155           continue;
156         }
157       }
158
159       /*
160       if (detectorIndex[0]<10) ftxt  << "  ->  " << detectorIndex[0];
161       else                     ftxt  << "  -> " << detectorIndex[0];
162       ftxt << "  " << detectorIndex[1];
163       if (detectorIndex[2]<10) ftxt << "  " << detectorIndex[2];
164       else                     ftxt << " " << detectorIndex[2];
165       ftxt << "  " << detectorIndex[3];
166       if (detectorIndex[4]<10) ftxt << "  " << detectorIndex[4];
167       else                     ftxt << " " << detectorIndex[4];
168
169       if (tofRawDatum->GetTOT()<10)                                            ftxt << "        " << tofRawDatum->GetTOT();
170       else if (tofRawDatum->GetTOT()>=10 && tofRawDatum->GetTOT()<100)         ftxt << "       " << tofRawDatum->GetTOT();
171       else if (tofRawDatum->GetTOT()>=100 && tofRawDatum->GetTOT()<1000)       ftxt << "      " << tofRawDatum->GetTOT();
172       else if (tofRawDatum->GetTOT()>=1000 && tofRawDatum->GetTOT()<10000)     ftxt << "     " << tofRawDatum->GetTOT();
173       else if (tofRawDatum->GetTOT()>=10000 && tofRawDatum->GetTOT()<100000)   ftxt << "    " << tofRawDatum->GetTOT();
174       else if (tofRawDatum->GetTOT()>=100000 && tofRawDatum->GetTOT()<1000000) ftxt << "   " << tofRawDatum->GetTOT();
175       else                                                                     ftxt << "  " << tofRawDatum->GetTOT();
176       if (tofRawDatum->GetTOF()<10)                                            ftxt << "        " << tofRawDatum->GetTOF() << endl;
177       else if (tofRawDatum->GetTOF()>=10 && tofRawDatum->GetTOF()<100)         ftxt << "       " << tofRawDatum->GetTOF() << endl;
178       else if (tofRawDatum->GetTOF()>=100 && tofRawDatum->GetTOF()<1000)       ftxt << "      " << tofRawDatum->GetTOF() << endl;
179       else if (tofRawDatum->GetTOF()>=1000 && tofRawDatum->GetTOF()<10000)     ftxt << "     " << tofRawDatum->GetTOF() << endl;
180       else if (tofRawDatum->GetTOF()>=10000 && tofRawDatum->GetTOF()<100000)   ftxt << "    " << tofRawDatum->GetTOF() << endl;
181       else if (tofRawDatum->GetTOF()>=100000 && tofRawDatum->GetTOF()<1000000) ftxt << "   " << tofRawDatum->GetTOF() << endl;
182       else                                                                     ftxt << "  " << tofRawDatum->GetTOF() << endl;
183       */
184
185       TClonesArray &aDigits = *tofDigits;
186       Int_t last = tofDigits->GetEntriesFast();
187
188       fTOFdigitMap->AddDigit(detectorIndex, last);
189
190       AliDebug(2,Form(" %3i -> %2i %2i %2i %2i %2i   %i  %i\n",
191                       last, detectorIndex[0], detectorIndex[1],
192                       detectorIndex[2], detectorIndex[4], detectorIndex[3],
193                       digit[1], digit[0]));
194
195       Int_t tracknum[3]={-1,-1,-1};
196       new (aDigits[last]) AliTOFdigit(tracknum, detectorIndex, digit);
197
198     } // while loop
199
200     clonesRawData->Clear();
201
202   } // DDL Loop
203
204   fTree->Fill();
205
206   //ftxt.close();
207
208   //delete cableLength;
209
210 }
211
212
213 /* ******************************************************* */
214 void AliEveTOFDigitsInfo::LoadDigits()
215 {
216   //
217   // Load TOF digits
218   //
219
220   TClonesArray *digitsTOF = 0x0;
221   AliTOFdigit *digs;
222
223   fTree->SetBranchAddress("TOF",&digitsTOF);
224   fTree->GetEntry(0);
225
226   Int_t vol[5] = {-1,-1,-1,-1,-1};
227
228   for (Int_t digitNumber=0; digitNumber<digitsTOF->GetEntries(); digitNumber++) {
229
230     if (digitNumber==digitsTOF->GetEntries()-1)
231       AliDebug(2,Form(" Hello  4 -> %3i digit of %i \n", digitNumber+1, digitsTOF->GetEntries()));
232
233     digs = (AliTOFdigit*)digitsTOF->UncheckedAt(digitNumber);
234
235     vol[0] = digs->GetSector(); // Sector Number (0-17)
236     vol[1] = digs->GetPlate();  // Plate Number (0-4)
237     vol[2] = digs->GetStrip();  // Strip Number (0-14/18)
238     vol[3] = digs->GetPadx();   // Pad Number in x direction (0-47)
239     vol[4] = digs->GetPadz();   // Pad Number in z direction (0-1)
240
241     fTOFdigitMap->AddDigit(vol, digitNumber);
242     if (digitNumber==digitsTOF->GetEntries()-1)
243       AliDebug(2,Form(" I am inside LoadDigits %3i \n", digitNumber));
244
245   }
246
247 }
248
249 /* ******************************************************* */
250
251 void AliEveTOFDigitsInfo::GetDigits(Int_t nSector, Int_t nPlate,
252                                     Int_t nStrip, Int_t nPadZ, Int_t nPadX,
253                                     Int_t indexDigit[3])
254 {
255   //
256   // Get TOF digit indices in the TOF volume
257   // (nSector, nPlate,nStrip,nPadZ,nPadX)
258   //
259
260   Int_t vol[5] = {nSector,nPlate,nStrip,nPadX,nPadZ};
261
262   fTOFdigitMap->GetDigitIndex(vol, indexDigit);
263   //for (Int_t ii=1; ii<3; ii++) indexDigit[ii]=-1;
264
265 }
266 /* ******************************************************* */
267
268 TClonesArray* AliEveTOFDigitsInfo::GetDigits(Int_t nSector, Int_t nPlate,
269                                              Int_t nStrip)
270 {
271   //
272   // Get TOF digits in the TOF volume
273   // (nSector, nPlate,nStrip)
274   //
275
276   Int_t newCounter = 0;
277   Int_t nDigitsInVolume[3] = {-1, -1, -1};
278   Int_t dummy[3] = {-1, -1, -1};
279   Int_t informations[4] = {-1, -1, -1, -1};
280
281   TClonesArray* digitsTOFnew = new TClonesArray("AliTOFdigit",  300);
282   TClonesArray &ldigits = *digitsTOFnew;
283
284   AliTOFdigit *digs;
285
286   TClonesArray *digitsTOF = 0x0;
287   fTree->SetBranchAddress("TOF",&digitsTOF);
288   fTree->GetEntry(0);
289
290   Int_t vol[5] = {nSector,nPlate,nStrip,-1,-1};
291
292   for(Int_t iPadZ=0; iPadZ<fGeom->NpadZ(); iPadZ++){
293     vol[4] = iPadZ;
294     for(Int_t iPadX=0; iPadX<fGeom->NpadX(); iPadX++) {
295       vol[3] = iPadX;
296
297       fTOFdigitMap->GetDigitIndex(vol, nDigitsInVolume);
298
299       for (Int_t ii=0; ii<3; ii++) {
300         //if (ii!=0) continue;
301         if (nDigitsInVolume[ii]>=0 ) {
302           AliDebug(2,Form("  nDigitsInVolume[%2i]  = %3i\n ", ii, nDigitsInVolume[ii]));
303           digs = (AliTOFdigit*)digitsTOF->UncheckedAt(nDigitsInVolume[ii]);
304           informations[0] = digs->GetTdc();
305           informations[1] = digs->GetAdc();
306           informations[2] = digs->GetToT();
307           informations[3] = digs->GetTdcND();
308           for(Int_t kk=0; kk<3; kk++) dummy[kk] = digs->GetTrack(kk);
309           new (ldigits[newCounter++]) AliTOFdigit(dummy, vol, informations);
310         }
311
312       }
313
314       for (Int_t ii=0; ii<4; ii++) informations[ii]=-1;
315       for (Int_t ii=0; ii<3; ii++) dummy[ii]=-1;
316       for (Int_t ii=0; ii<3; ii++) nDigitsInVolume[ii]=-1;
317
318     }
319   }
320
321   if (digitsTOFnew)
322     AliDebug(2, Form("Sector %2i   Plate %1i  Strip %2i  -> number of digits %3i \n",
323                      nSector, nPlate, nStrip, digitsTOFnew->GetEntries()));
324
325   return digitsTOFnew;
326
327 }
328 /* ******************************************************* */
329
330 TClonesArray* AliEveTOFDigitsInfo::GetDigits(Int_t nSector)
331 {
332   //
333   // Get TOF digits in the TOF SM nSector
334   //
335
336   const Int_t kND = AliTOFDigitMap::kMaxDigitsPerPad;
337
338   Int_t newCounter = 0;
339   Int_t nDigitsInVolume[kND];
340   Int_t dummy[3];
341   Int_t informations[4];
342
343   Int_t nStrips=19;
344
345   TClonesArray* digitsTOFnew = new TClonesArray("AliTOFdigit",  300);
346   TClonesArray &ldigits = *digitsTOFnew;
347
348   AliTOFdigit *digs;
349
350   TClonesArray *digitsTOF = 0x0;
351   fTree->SetBranchAddress("TOF",&digitsTOF);
352   fTree->GetEntry(0);
353
354   //Int_t nSector = 1;
355   Int_t vol[5] = {nSector,-1,-1,-1,-1};
356  
357   for(Int_t iPlate=0; iPlate<fGeom->NPlates(); iPlate++){
358     vol[1] = iPlate;
359     if(iPlate==2) nStrips=15;
360     else nStrips=19;
361
362     for(Int_t iStrip=0; iStrip<nStrips; iStrip++){
363       vol[2] = iStrip;
364
365       for(Int_t iPadZ=0; iPadZ<fGeom->NpadZ(); iPadZ++){
366         vol[4] = iPadZ;
367
368         for(Int_t iPadX=0; iPadX<fGeom->NpadX(); iPadX++) {
369
370           for (Int_t ii=0; ii<4; ii++) informations[ii]=-1;
371           for (Int_t ii=0; ii<3; ii++) dummy[ii]=-1;
372           for (Int_t ii=0; ii<kND; ii++) nDigitsInVolume[ii]=-1;
373
374           vol[3] = iPadX;
375
376           fTOFdigitMap->GetDigitIndex(vol, nDigitsInVolume);
377
378           for (Int_t ii=0; ii<kND; ii++) {
379             //if (ii!=0) continue;
380             if (nDigitsInVolume[ii]>=0 ) {
381
382               digs = (AliTOFdigit*)digitsTOF->UncheckedAt(nDigitsInVolume[ii]);
383               informations[0] = digs->GetTdc();
384               informations[1] = digs->GetAdc();
385               informations[2] = digs->GetToT();
386               informations[3] = digs->GetTdcND();
387               for(Int_t kk=0; kk<3; kk++) dummy[kk] = digs->GetTrack(kk);
388               new (ldigits[newCounter++]) AliTOFdigit(dummy, vol, informations);
389
390               AliDebug(2,Form(" %2i -> %2i %2i %2i %2i %2i %7i %7i\n",
391                               nDigitsInVolume[ii],
392                               vol[0], vol[1], vol[2], vol[4], vol[3],
393                               informations[1], informations[0]));
394
395             }
396
397           }
398
399         }
400       }
401     }
402   }
403
404   if (digitsTOFnew)
405     AliDebug(2,Form("Sector %2i  -> number of digits %3i \n",
406                     nSector, digitsTOFnew->GetEntries()));
407
408   return digitsTOFnew;
409
410 }
411 /* ******************************************************* */
412
413 Int_t AliEveTOFDigitsInfo::GetTOFInfos() const
414 {
415   //
416   // Return number of TOF digits
417   //
418
419   return fTOFdigitMap->GetFilledCellNumber();
420
421 }
422
423 /* ******************************************************* */
424 Int_t AliEveTOFDigitsInfo::IsStripFilled(Int_t iSector, Int_t iPlate, Int_t iStrip)
425 {
426   //
427   // Return number of TOF digits
428   // in volume (iSector,iPlate,iStrip)
429   //
430
431   Int_t vol[5] = {iSector, iPlate, iStrip, -1, -1};
432
433   Int_t index = 0;
434
435   for (Int_t iPadZ=0; iPadZ<fGeom->NpadZ(); iPadZ++) 
436     for (Int_t iPadX=0; iPadX<fGeom->NpadX(); iPadX++) 
437       {
438         vol[3] = iPadX;
439         vol[4] = iPadZ;
440         if (fTOFdigitMap->GetDigitIndex(vol,0)>=0) index++;
441       }
442
443   return index;
444
445 }
446
447 /* ******************************************************* */
448 /*
449 void AliEveTOFDigitsInfo::GetDigits()
450 {
451
452   for (Int_t iSector=0; iSector<fGeom->NSectors(); iSector++) {
453
454     fNewTree = new TTree();
455
456
457
458
459   }
460
461 }
462 */