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