]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitsManager.cxx
7c759249e7db8bbef5cfb516de8278ae3361ef5e
[u/mrichter/AliRoot.git] / TRD / AliTRDdigitsManager.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  Manages the digits and the track dictionary in the form of               //
21 //  AliTRDdataArray objects.                                                 //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include <Riostream.h>
26  
27 #include <TROOT.h>
28 #include <TTree.h>                                                              
29 #include <TFile.h>
30
31 #include "AliRun.h"
32 #include "AliLog.h"
33
34 #include "AliTRDdigitsManager.h"
35 #include "AliTRDsegmentArray.h"
36 #include "AliTRDdataArrayI.h"
37 #include "AliTRDdigit.h"
38 #include "AliTRDgeometry.h"
39
40 #include "AliTRDSignalIndex.h"
41
42 ClassImp(AliTRDdigitsManager)
43
44 //_____________________________________________________________________________
45
46   // Number of track dictionary arrays
47   const Int_t AliTRDdigitsManager::fgkNDict = kNDict;
48
49 //_____________________________________________________________________________
50 AliTRDdigitsManager::AliTRDdigitsManager()
51   :TObject()
52   ,fEvent(0)
53   ,fTree(0)
54   ,fDigits(0)
55   ,fIsRaw(0)
56   ,fSDigits(0)
57   ,fSignalIndexes(NULL)
58 {
59   //
60   // Default constructor
61   //
62
63   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
64     fDictionary[iDict] = NULL;
65   }
66   
67   //fSignalIndexes = new TList();
68   fSignalIndexes = new TObjArray(AliTRDgeometry::Ndet());
69   
70 }
71
72 //_____________________________________________________________________________
73 AliTRDdigitsManager::AliTRDdigitsManager(const AliTRDdigitsManager &m)
74   :TObject(m)
75   ,fEvent(m.fEvent)
76   ,fTree(0)
77   ,fDigits(0)
78   ,fIsRaw(m.fIsRaw)
79   ,fSDigits(m.fSDigits)
80   ,fSignalIndexes(NULL)
81 {
82   //
83   // AliTRDdigitsManager copy constructor
84   //
85
86 }
87
88 //_____________________________________________________________________________
89 AliTRDdigitsManager::~AliTRDdigitsManager()
90 {
91   //
92   // AliTRDdigitsManager destructor
93   //
94
95   if (fDigits) {
96     fDigits->Delete();
97     delete fDigits;
98     fDigits = NULL;
99   }
100
101   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
102     fDictionary[iDict]->Delete();
103     delete fDictionary[iDict];
104     fDictionary[iDict] = NULL;
105   }
106
107   delete fSignalIndexes;
108   fSignalIndexes = NULL;
109 //   for (Int_t i = 0; i < AliTRDgeometry::Ndet(); i++)
110 //     delete fSignalIndexes[i];
111
112 }
113
114 //_____________________________________________________________________________
115 AliTRDdigitsManager &AliTRDdigitsManager::operator=(const AliTRDdigitsManager &m)
116 {
117   //
118   // Assignment operator
119   //
120
121   if (this != &m) ((AliTRDdigitsManager &) m).Copy(*this);
122   return *this;
123
124 }
125
126 //_____________________________________________________________________________
127 void AliTRDdigitsManager::Copy(TObject &m) const
128 {
129   //
130   // Copy function
131   //
132
133   ((AliTRDdigitsManager &) m).fIsRaw   = fIsRaw;
134   ((AliTRDdigitsManager &) m).fEvent   = fEvent;
135   ((AliTRDdigitsManager &) m).fSDigits = fSDigits;
136   
137   ((AliTRDdigitsManager &) m).fSignalIndexes = fSignalIndexes;
138   TObject::Copy(m);
139
140 }
141
142 //_____________________________________________________________________________
143 void AliTRDdigitsManager::CreateArrays()
144 {
145   //
146   // Create the data arrays
147   //
148
149   fDigits = new AliTRDsegmentArray("AliTRDdataArrayI",AliTRDgeometry::Ndet());
150
151   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
152     fDictionary[iDict] = new AliTRDsegmentArray("AliTRDdataArrayI"
153                                                ,AliTRDgeometry::Ndet());
154   }
155
156   for (Int_t i = 0; i < AliTRDgeometry::Ndet(); i++)
157     {
158       fSignalIndexes->AddLast(new AliTRDSignalIndex());
159     }
160 }
161 //_____________________________________________________________________________
162 void AliTRDdigitsManager::ResetArrays()
163 {
164   //
165   // Reset the data arrays
166   //
167
168   if (fDigits) {
169     delete fDigits;
170   }
171   fDigits = new AliTRDsegmentArray("AliTRDdataArrayI",AliTRDgeometry::Ndet());
172
173   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
174     if (fDictionary[iDict]) {  
175       delete fDictionary[iDict];
176     }
177     fDictionary[iDict] = new AliTRDsegmentArray("AliTRDdataArrayI"
178                                                ,AliTRDgeometry::Ndet());
179   }
180
181   for (Int_t i = 0; i < AliTRDgeometry::Ndet(); i++)
182     {
183       AliTRDSignalIndex *idx = (AliTRDSignalIndex *)fSignalIndexes->At(i);
184       idx->Reset();
185     }
186 }
187
188 //_____________________________________________________________________________
189 void AliTRDdigitsManager::SetRaw()
190 {
191   //
192   // Switch on the raw digits flag
193   //
194
195   fIsRaw = kTRUE;
196   if (fDigits)
197     fDigits->SetBit(AliTRDdigit::RawDigit());
198   
199 }
200
201 //_____________________________________________________________________________
202 Short_t AliTRDdigitsManager::GetDigitAmp(Int_t row, Int_t col,Int_t time
203                                        , Int_t det) const
204 {
205   //
206   // Returns the amplitude of a digit
207   //
208
209   if (!GetDigits(det)) return 0;
210   return ((Short_t) GetDigits(det)->GetData(row,col,time));
211
212 }
213  
214 //_____________________________________________________________________________
215 Bool_t AliTRDdigitsManager::MakeBranch(TTree *tree)
216 {
217   //
218   // Creates the tree and branches for the digits and the dictionary
219   //
220
221   Int_t buffersize = 64000;
222
223   Bool_t status = kTRUE;
224
225   if (tree) {
226     fTree = tree;
227   }
228
229   // Make the branch for the digits
230   if (fDigits) {
231     const AliTRDdataArray *kDigits = (AliTRDdataArray *) fDigits->At(0);
232     if (kDigits) {
233       if (!fTree) return kFALSE;
234       TBranch* branch = fTree->GetBranch("TRDdigits");
235       if (!branch) fTree->Branch("TRDdigits",kDigits->IsA()->GetName(),
236                                  &kDigits,buffersize,99);
237       AliDebug(1,"Making branch TRDdigits\n");
238     }
239     else {
240       status = kFALSE;
241     }
242   }
243   else {
244     status = kFALSE;
245   }
246
247   // Make the branches for the dictionaries
248   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
249     Char_t branchname[15];
250     sprintf(branchname,"TRDdictionary%d",iDict);
251     if (fDictionary[iDict]) {
252       const AliTRDdataArray *kDictionary = 
253               (AliTRDdataArray *) fDictionary[iDict]->At(0);
254       if (kDictionary) {
255         if (!fTree) return kFALSE;
256         TBranch* branch = fTree->GetBranch(branchname);
257         if (!branch) fTree->Branch(branchname,kDictionary->IsA()->GetName(),
258                                    &kDictionary,buffersize,99);
259         AliDebug(1,Form("Making branch %s\n",branchname));
260       }
261       else {
262         status = kFALSE;
263       }
264     }
265     else {
266       status = kFALSE;
267     }
268   }
269
270   return status;
271
272 }
273
274 //_____________________________________________________________________________
275 Bool_t AliTRDdigitsManager::ReadDigits(TTree *tree)
276 {
277   //
278   // Reads the digit information from the input file
279   //
280
281   Bool_t status = kTRUE;
282
283   if (tree) {
284
285     fTree = tree;
286
287   }
288
289   if (!fDigits) {
290     AliDebug(1,"Create the data arrays.\n");
291     CreateArrays();
292   }
293
294   status = fDigits->LoadArray("TRDdigits",fTree);
295
296   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
297     Char_t branchname[15];
298     sprintf(branchname,"TRDdictionary%d",iDict);
299     status = fDictionary[iDict]->LoadArray(branchname,fTree);
300   }  
301
302   if (fDigits->TestBit(AliTRDdigit::RawDigit())) {
303     fIsRaw = kTRUE;
304   }
305   else {
306     fIsRaw = kFALSE;
307   }
308
309   return kTRUE;
310
311 }
312
313 //_____________________________________________________________________________
314 Bool_t AliTRDdigitsManager::WriteDigits()
315 {
316   //
317   // Writes out the TRD-digits and the dictionaries
318   //
319
320   // Store the contents of the segment array in the tree
321   if (!fDigits->StoreArray("TRDdigits",fTree)) {
322     AliError("Error while storing digits in branch TRDdigits\n");
323     return kFALSE;
324   }
325   for (Int_t iDict = 0; iDict < kNDict; iDict++) {
326     Char_t branchname[15];
327     sprintf(branchname,"TRDdictionary%d",iDict);
328     if (!fDictionary[iDict]->StoreArray(branchname,fTree)) {
329       AliError(Form("Error while storing dictionary in branch %s\n",branchname));
330       return kFALSE;
331     }
332   }
333
334   // Write the new tree to the output file
335   fTree->AutoSave();  // Modification by Jiri
336
337   return kTRUE;
338
339 }
340
341 //_____________________________________________________________________________
342 AliTRDdigit *AliTRDdigitsManager::GetDigit(Int_t row, Int_t col
343                                          , Int_t time, Int_t det) const
344 {
345   // 
346   // Creates a single digit object 
347   //
348
349   Int_t digits[4];
350   Int_t amp[1];
351
352   digits[0] = det;
353   digits[1] = row;
354   digits[2] = col;
355   digits[3] = time;
356
357   amp[0]    = GetDigits(det)->GetData(row,col,time);
358   
359   return (new AliTRDdigit(fIsRaw,digits,amp));
360
361 }
362
363 //_____________________________________________________________________________
364 Int_t AliTRDdigitsManager::GetTrack(Int_t track
365                                   , Int_t row, Int_t col, Int_t time
366                                   , Int_t det) const
367 {
368   // 
369   // Returns the MC-track numbers from the dictionary.
370   //
371
372   if ((track < 0) || (track >= kNDict)) {
373     AliError(Form("track %d out of bounds (size: %d, this: 0x%08x)"
374                  ,track,kNDict,this));
375     return -1;
376   }
377
378   // Array contains index+1 to allow data compression
379   return (GetDictionary(det,track)->GetData(row,col,time) - 1);
380
381 }
382
383 //_____________________________________________________________________________
384 AliTRDdataArrayI *AliTRDdigitsManager::GetDigits(Int_t det) const
385 {
386   //
387   // Returns the digits array for one detector
388   //
389
390   if (!fDigits) return 0x0;
391   return (AliTRDdataArrayI *) fDigits->At(det);
392
393 }
394
395 //_____________________________________________________________________________
396 AliTRDdataArrayI *AliTRDdigitsManager::GetDictionary(Int_t det, Int_t i) const
397 {
398   //
399   // Returns the dictionary for one detector
400   //
401
402   return (AliTRDdataArrayI *) fDictionary[i]->At(det);
403
404 }
405
406 //_____________________________________________________________________________
407 Int_t AliTRDdigitsManager::GetTrack(Int_t track, AliTRDdigit *Digit) const
408 {
409   // 
410   // Returns the MC-track numbers from the dictionary for a given digit
411   //
412
413   Int_t row  = Digit->GetRow();
414   Int_t col  = Digit->GetCol();
415   Int_t time = Digit->GetTime();
416   Int_t det  = Digit->GetDetector();
417
418   return GetTrack(track,row,col,time,det);
419
420 }
421
422 //_____________________________________________________________________________
423 AliTRDSignalIndex *AliTRDdigitsManager::GetIndexes(Int_t det) 
424 {
425   // 
426   // Returns indexes of active pads
427   //
428
429   return (AliTRDSignalIndex*)fSignalIndexes->At(det);
430
431 }
432
433 //_____________________________________________________________________________
434 void AliTRDdigitsManager::RemoveDigits(Int_t det) 
435 {
436   // 
437   // Clear memory
438   //
439
440   fDigits->ClearSegment(det);
441
442 }
443
444 //_____________________________________________________________________________
445 void AliTRDdigitsManager::RemoveDictionaries(Int_t det) 
446 {
447   // 
448   // Clear memory
449   //
450
451   for (Int_t i = 0; i < kNDict; i++) {
452     fDictionary[i]->ClearSegment(det);
453   }
454
455 }
456
457 //_____________________________________________________________________________
458 void AliTRDdigitsManager::ClearIndexes(Int_t det) 
459 {
460   // 
461   // Clear memory
462   //
463   fSignalIndexes->At(det)->Clear();  
464 }
465
466
467 //_____________________________________________________________________________
468 Bool_t AliTRDdigitsManager::BuildIndexes(Int_t det)
469 {
470   //
471   // Build the list of indices
472   //
473
474   Int_t nRows = 0;
475   Int_t nCols = 0;
476   Int_t nTbins = 0;
477
478   AliTRDgeometry    geom;
479   AliTRDdataArrayI *digits = GetDigits(det);
480   //digits should be expanded by now!!!
481   if (digits->GetNtime() > 0) 
482     {
483       digits->Expand();
484       nRows = digits->GetNrow();
485       nCols = digits->GetNcol();
486       nTbins = digits->GetNtime();
487
488       //AliInfo(Form("rows %d cols %d tbins %d", nRows, nCols, nTbins));
489
490       AliTRDSignalIndex* indexes = GetIndexes(det);
491       indexes->SetSM(geom.GetSector(det));
492       indexes->SetChamber(geom.GetChamber(det));
493       indexes->SetPlane(geom.GetPlane(det));
494       indexes->SetDetNumber(det);
495       if (indexes->IsAllocated() == kFALSE)
496         {
497           indexes->Allocate(nRows, nCols, nTbins);
498           //AliInfo(Form("Allocating 0x%x %d", indexes->GetArray(), indexes->GetArray()->GetSize()));
499         }
500       for (Int_t ir = 0; ir < nRows; ir++)
501         {
502           for (Int_t ic = 0; ic < nCols; ic++)
503             {
504               for (Int_t it = 0; it < nTbins; it++)
505                 {         
506                   //AliInfo(Form("row %d col %d tbin %d", ir, ic, it));
507                   
508                   Int_t isig = digits->GetDataUnchecked(ir, ic, it);
509                   if (isig > 0)
510                     {
511                       //AliInfo(Form("row %d col %d tbin %d", ir, ic, it));
512                       indexes->AddIndexTBin(ir, ic, it);            
513                     }
514                 } //tbins
515             } //cols
516         } // rows
517     } // if GetNtime
518   else
519     {
520       return kFALSE;
521     }
522   return kTRUE;
523
524 }