]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitsManager.cxx
1d36aedb243e8197b36d77d572b2934b2b5e7de7
[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 //  TObjArray objects                                                        //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include <TTree.h>                                                              
26
27 #include "AliLog.h"
28
29 #include "AliTRDdigitsManager.h"
30 #include "AliTRDarrayDictionary.h"
31 #include "AliTRDarrayADC.h"
32 #include "AliTRDarraySignal.h"
33 #include "AliTRDdigit.h"
34 #include "AliTRDdigitsParam.h"
35 #include "AliTRDSimParam.h"
36 #include "AliTRDgeometry.h"
37 #include "AliTRDSignalIndex.h"
38
39 ClassImp(AliTRDdigitsManager)
40
41 //_____________________________________________________________________________
42
43   // Number of track dictionary arrays
44   const Int_t AliTRDdigitsManager::fgkNDict = kNDict;
45
46 //_____________________________________________________________________________
47 AliTRDdigitsManager::AliTRDdigitsManager(Bool_t rawRec)
48   :TObject()
49   ,fEvent(0)
50   ,fTree(0)
51   ,fDigits(0) 
52   ,fHasSDigits(0)
53   ,fSignalIndexes(NULL)
54   ,fUseDictionaries(kTRUE)
55   ,fDets(AliTRDgeometry::Ndet())
56   ,fRawRec(rawRec)
57   ,fDigitsParam(0)
58 {
59   //
60   // Default constructor
61   //
62   
63   if (fRawRec)
64     {
65       fDets   = 1;
66       fRawRec = kTRUE;
67     }
68
69   for (Int_t iDict = 0; iDict < kNDict; iDict++) 
70     {
71       fDict[iDict] = NULL;
72     }
73   
74 }
75
76 //_____________________________________________________________________________
77 AliTRDdigitsManager::AliTRDdigitsManager(const AliTRDdigitsManager &m)
78   :TObject(m)
79   ,fEvent(m.fEvent)
80   ,fTree(0)
81   ,fDigits(0) 
82   ,fHasSDigits(m.fHasSDigits)
83   ,fSignalIndexes(NULL)
84   ,fUseDictionaries(kTRUE)
85   ,fDets(m.fDets)
86   ,fRawRec(m.fRawRec)
87   ,fDigitsParam(NULL)
88 {
89   //
90   // AliTRDdigitsManager copy constructor
91   //
92
93 }
94
95 //_____________________________________________________________________________
96 AliTRDdigitsManager::~AliTRDdigitsManager()
97 {
98   //
99   // AliTRDdigitsManager destructor
100   //
101
102
103   if (fDigits) 
104     {
105       fDigits->Delete();
106       delete fDigits;
107       fDigits = NULL;
108     }
109
110   for (Int_t iDict = 0; iDict < kNDict; iDict++) 
111     {
112       if(fDict[iDict])
113         {
114           fDict[iDict]->Delete();
115           delete fDict[iDict];
116           fDict[iDict] = NULL;
117         }
118     }
119
120   if (fSignalIndexes) 
121     {
122       fSignalIndexes->Delete();
123       delete fSignalIndexes;
124       fSignalIndexes = NULL;
125     }
126
127   if (fDigitsParam)
128     {
129       delete fDigitsParam;
130       fDigitsParam = NULL;
131     }
132
133 }
134
135 //_____________________________________________________________________________
136 AliTRDdigitsManager &AliTRDdigitsManager::operator=(const AliTRDdigitsManager &m)
137 {
138   //
139   // Assignment operator
140   //
141
142   if (this != &m) 
143     {
144       ((AliTRDdigitsManager &) m).Copy(*this);
145     }
146
147   return *this;
148
149 }
150
151 //_____________________________________________________________________________
152 void AliTRDdigitsManager::Copy(TObject &m) const
153 {
154   //
155   // Copy function
156   //
157
158   ((AliTRDdigitsManager &) m).fEvent           = fEvent;
159   ((AliTRDdigitsManager &) m).fHasSDigits      = fHasSDigits;
160   ((AliTRDdigitsManager &) m).fDigits          = NULL;
161   for(Int_t i = 0; i < kNDict; i++)
162     {
163       ((AliTRDdigitsManager &) m).fDict[i]  = NULL;
164     }
165   ((AliTRDdigitsManager &) m).fSignalIndexes   = NULL;
166   ((AliTRDdigitsManager &) m).fUseDictionaries = fUseDictionaries;
167   ((AliTRDdigitsManager &) m).fDets            = fDets;
168   ((AliTRDdigitsManager &) m).fRawRec          = fRawRec;
169   ((AliTRDdigitsManager &) m).fDigitsParam     = NULL;
170
171   TObject::Copy(m);
172
173 }
174
175 //_____________________________________________________________________________
176 void AliTRDdigitsManager::CreateArrays()
177 {
178   //
179   // Create the data arrays
180   //
181
182   if (fHasSDigits) 
183     {
184       if (fDigits)                                        
185         {                                                   
186           fDigits->Delete();                                
187           delete fDigits;                                   
188         }                                                    
189       fDigits = new TObjArray(fDets);
190       for (Int_t index = 0; index < fDets; index++) 
191         {
192           fDigits->AddAt(new AliTRDarraySignal(),index);
193         }
194     }
195   else 
196     {
197       if (fDigits)                                          
198         {                                                    
199           fDigits->Delete();                                
200           delete fDigits;                                   
201         }                                                   
202       fDigits = new TObjArray(fDets);    
203       for (Int_t index = 0; index < fDets; index++) 
204         {
205           fDigits->AddAt(new AliTRDarrayADC(),index);
206         }
207     }
208
209   if (fUseDictionaries) 
210     {
211       for (Int_t iDict = 0; iDict < kNDict; iDict++)
212         {
213           if (fDict[iDict])                                   
214             {
215               fDict[iDict]->Delete();                                
216               delete fDict[iDict];                                    
217             }
218           fDict[iDict] = new TObjArray(fDets);
219           for (Int_t index = 0; index < fDets; index++) 
220             {
221               fDict[iDict]->AddAt(new AliTRDarrayDictionary(),index);
222             }
223         }
224     }
225   
226   if (fSignalIndexes)
227     {
228       fSignalIndexes->Delete();
229       delete fSignalIndexes;
230     }
231   fSignalIndexes = new TObjArray(fDets);
232   for (Int_t i = 0; i < fDets; i++)
233     {
234       fSignalIndexes->AddLast(new AliTRDSignalIndex());
235     }
236
237   if (fDigitsParam)
238     {
239       delete fDigitsParam;
240     }
241   fDigitsParam = new AliTRDdigitsParam();
242   fDigitsParam->SetNTimeBins(AliTRDSimParam::Instance()->GetNTimeBins());
243   fDigitsParam->SetADCbaseline(AliTRDSimParam::Instance()->GetADCbaseline());
244
245 }
246
247 //_____________________________________________________________________________
248 void AliTRDdigitsManager::ResetArrays()
249 {
250   //
251   // Reset the data arrays
252   //
253
254   if (fDigits)
255     {
256       fDigits->Delete();
257       delete fDigits;
258     }
259   if (fHasSDigits)
260     {
261       fDigits = new TObjArray(fDets);     
262       for (Int_t index = 0; index < fDets; index++) 
263         {
264           fDigits->AddAt(new AliTRDarraySignal(),index);
265         }
266     }
267   else
268     {
269       fDigits = new TObjArray(fDets);
270       for (Int_t index = 0; index < fDets; index++)
271         {
272           fDigits->AddAt(new AliTRDarrayADC(),index);
273         }
274     }
275   
276   for (Int_t iDict = 0; iDict < kNDict; iDict++)
277     {
278       if (fDict[iDict])
279         {
280           fDict[iDict]->Delete();
281           delete fDict[iDict];
282           fDict[iDict] = NULL;
283         }
284     }
285   if (fUseDictionaries) 
286     {
287       for (Int_t iDict = 0; iDict < kNDict; iDict++)
288         {
289           fDict[iDict] = new TObjArray(fDets);
290           for (Int_t index = 0; index < fDets; index++)
291             {
292               fDict[iDict]->AddAt(new AliTRDarrayDictionary(),index);
293             }
294         }
295     }
296   
297   if (fSignalIndexes)
298     {
299       fSignalIndexes->Delete();
300       delete fSignalIndexes;
301     }
302   fSignalIndexes = new TObjArray(fDets);
303   for (Int_t i = 0; i < fDets; i++)
304     {
305       fSignalIndexes->AddLast(new AliTRDSignalIndex());
306     }
307
308 }
309
310 //_____________________________________________________________________________
311 void AliTRDdigitsManager::ResetArrays(Int_t det)
312 {
313   //
314   // Reset the data arrays
315   //
316
317   Int_t recoDet = fRawRec ? 0 : det;
318
319   RemoveDigits(recoDet);
320   RemoveDictionaries(recoDet);
321   RemoveIndexes(recoDet);
322
323   if (fHasSDigits)
324     {
325       fDigits->AddAt(new AliTRDarraySignal(),recoDet);
326     }
327   else
328     {
329       fDigits->AddAt(new AliTRDarrayADC(),recoDet);
330     }
331
332   if (fUseDictionaries) 
333     {
334       for (Int_t iDict = 0; iDict < kNDict; iDict++)
335         {
336           fDict[iDict]->AddAt(new AliTRDarrayDictionary(),recoDet);
337         }
338     }
339   
340   fSignalIndexes->AddAt(new AliTRDSignalIndex(),recoDet);
341
342 }
343
344 //_____________________________________________________________________________
345 void AliTRDdigitsManager::ClearArrays(Int_t det)
346 {
347   //
348   // Reset the data arrays
349   //
350
351   Int_t recoDet = fRawRec ? 0 : det;
352
353   if (fHasSDigits)
354     {
355       ((AliTRDarraySignal*)fDigits->At(recoDet))->Reset();
356     }
357   else
358     {
359       ((AliTRDarrayADC*)fDigits->At(recoDet))->ConditionalReset((AliTRDSignalIndex*)fSignalIndexes->At(recoDet));
360     }
361
362   if (fUseDictionaries) 
363     {
364       for (Int_t iDict = 0; iDict < kNDict; iDict++)
365         {
366           ((AliTRDarrayDictionary*)fDict[iDict]->At(recoDet))->Reset();
367         }
368     }
369   
370   ((AliTRDSignalIndex*)fSignalIndexes->At(recoDet))->ResetContent();
371
372 }
373
374 //_____________________________________________________________________________
375 Short_t AliTRDdigitsManager::GetDigitAmp(Int_t row, Int_t col,Int_t time, Int_t det) const
376 {
377   //
378   // Returns the amplitude of a digit
379   //
380
381   if (!GetDigits(det)) return 0;
382   
383   return ((Short_t) ((AliTRDarrayADC *) GetDigits(det))->GetDataBits(row,col,time));
384
385 }
386
387 //_____________________________________________________________________________
388 UChar_t AliTRDdigitsManager::GetPadStatus(Int_t row, Int_t col, Int_t time, Int_t det) const
389 {
390   //
391   // Returns the pad status for the requested pad
392   //
393         
394   if (!GetDigits(det)) return 0;
395
396   return ((UChar_t) ((AliTRDarrayADC *) GetDigits(det))->GetPadStatus(row,col,time));
397  
398 }
399
400 //_____________________________________________________________________________
401 Bool_t AliTRDdigitsManager::MakeBranch(TTree * const tree)  
402 {
403   //
404   // Creates the tree and branches for the digits and the dictionary
405   //
406
407   Int_t  buffersize = 64000;
408   Bool_t status     = kTRUE;
409
410   if (tree) 
411     {
412       fTree = tree;
413     }
414
415   // Make the branch for the digits
416   if (fDigits) 
417     {
418       if (fHasSDigits)
419         {
420           const AliTRDarraySignal *kDigits = (AliTRDarraySignal *) fDigits->At(0); 
421           if (kDigits) 
422             {
423               if (!fTree) return kFALSE;
424               AliDebug(1,"Making branch for SDigits!\n");
425               TBranch *branch = fTree->GetBranch("TRDdigits");
426               if (!branch)
427                 {
428                   fTree->Branch("TRDdigits","AliTRDarraySignal",&kDigits,buffersize,99);
429                 }
430               AliDebug(1,"Making branch TRDdigits\n");
431             }
432           else 
433             {
434               status = kFALSE;
435             }
436         }
437
438       if (!fHasSDigits)
439         {
440           const AliTRDarrayADC *kDigits = (AliTRDarrayADC *) fDigits->At(0);
441           if (kDigits) 
442             {
443               if (!fTree) return kFALSE;
444               AliDebug(1,"Making branch for Digits!\n");
445               TBranch *branch = fTree->GetBranch("TRDdigits");
446               if (!branch) 
447                 {
448                   fTree->Branch("TRDdigits","AliTRDarrayADC",&kDigits,buffersize,99);
449                 }
450               AliDebug(1,"Making branch TRDdigits\n");        
451             }
452           else 
453             {
454               status = kFALSE;
455             }
456         }
457
458     }
459   else
460     {
461
462       status = kFALSE;
463
464     }
465   
466   if (fUseDictionaries) 
467     {
468       // Make the branches for the dictionaries
469       for (Int_t iDict = 0; iDict < kNDict; iDict++) 
470         {
471           Char_t branchname[15];
472           sprintf(branchname,"TRDdictionary%d",iDict); 
473           if (fDict[iDict]) 
474             {
475               const AliTRDarrayDictionary *kDictionary = (AliTRDarrayDictionary *) fDict[iDict]->At(0);
476               if (kDictionary) 
477                 {
478                   if (!fTree) return kFALSE;
479                   AliDebug(2,"Making branch for dictionary!\n");
480                   TBranch *branch = fTree->GetBranch(branchname);
481                   if (!branch) 
482                     {
483                       fTree->Branch(branchname,"AliTRDarrayDictionary",&kDictionary,buffersize,99);
484                     }
485                   AliDebug(1,Form("Making branch %s\n",branchname));
486                 }
487               else 
488                 {
489                   status = kFALSE;
490                 }
491             }
492           else 
493             {
494               status = kFALSE;
495             }
496         }
497     }
498
499   if (fDigitsParam)
500     {
501       const AliTRDdigitsParam *kDigitsParam = fDigitsParam;
502       if (!fTree) return kFALSE;
503       TBranch *branch = fTree->GetBranch("TRDdigitsParam");
504       if (!branch) 
505         {
506           fTree->Branch("TRDdigitsParam","AliTRDdigitsParam",&kDigitsParam,buffersize,99);
507         }
508       AliDebug(1,"Making branch AliTRDdigitsParam\n");
509     }
510
511   return status;
512
513 }
514
515 //_____________________________________________________________________________
516 Bool_t AliTRDdigitsManager::ReadDigits(TTree * const tree)
517 {
518   //
519   // Reads the digit information from the input file
520   //
521
522   Bool_t status = kTRUE;
523
524   if (tree) 
525     {
526       fTree = tree;
527     }
528
529   if (!fDigits) 
530     {
531       AliDebug(1,"Create the data arrays.\n");
532       CreateArrays();
533     }
534
535   status = LoadArrayDigits();
536
537   if (fUseDictionaries) 
538     {
539       status = LoadArrayDict();
540       if (status == kFALSE) 
541         {
542           fUseDictionaries = kFALSE;
543           AliWarning("Unable to load dict arrays. Will not use them.\n");
544         }
545     }  
546
547   if (!LoadDigitsParam()) {
548     AliWarning("Could not read digits parameter.");
549     if (fDigitsParam) {
550       delete fDigitsParam;
551     }
552     AliWarning(Form("Create default version of digits parameter (NTimeBin=%d).\n"
553                    ,AliTRDSimParam::Instance()->GetNTimeBins()));
554     fDigitsParam = new AliTRDdigitsParam();
555     fDigitsParam->SetNTimeBins(AliTRDSimParam::Instance()->GetNTimeBins());
556     fDigitsParam->SetADCbaseline(AliTRDSimParam::Instance()->GetADCbaseline());
557   }
558
559   return status;
560
561 }
562
563 //_____________________________________________________________________________
564 Bool_t AliTRDdigitsManager::WriteDigits()
565 {
566   //
567   // Writes out the TRD-digits and the dictionaries
568   //
569
570   if (!StoreArrayDigits())
571     {
572       AliError("Error while storing digits\n");
573       return kFALSE;
574     }
575
576   if (fUseDictionaries) 
577     {
578       if (!StoreArrayDict())
579         {
580           AliError("Error while storing dictionaries in branch\n");
581           return kFALSE;
582         }
583     }
584   
585   // Write the new tree to the output file
586   fTree->AutoSave();
587
588   return kTRUE;
589
590 }
591
592 //_____________________________________________________________________________
593 AliTRDdigit *AliTRDdigitsManager::GetDigit(Int_t row
594                                          , Int_t col
595                                          , Int_t time
596                                          , Int_t det) const
597 {
598   // 
599   // Creates a single digit object 
600   //
601
602   Int_t digits[4]; 
603   Int_t amp[1];
604
605   digits[0] = det;
606   digits[1] = row;
607   digits[2] = col;
608   digits[3] = time;
609
610   amp[0]    = ((AliTRDarrayADC *) GetDigits(det))->GetData(row,col,time);
611   
612   return (new AliTRDdigit(digits,amp));
613
614 }
615
616 //_____________________________________________________________________________
617 Int_t AliTRDdigitsManager::GetTrack(Int_t track
618                                   , Int_t row
619                                   , Int_t col
620                                   , Int_t time
621                                   , Int_t det) const
622 {
623   // 
624   // Returns the MC-track numbers from the dictionary.
625   //
626
627   if ((track < 0) || (track >= kNDict)) 
628     {
629       AliError(Form("track %d out of bounds (size: %d, this: 0x%08x)",track,kNDict,this));
630       return -1;
631     }
632
633   if (fUseDictionaries == kFALSE) 
634     {
635       return -1;
636     }
637
638   // Array contains index+1 to allow data compression--->Changed
639   return (((AliTRDarrayDictionary *) GetDictionary(det,track))->GetData(row,col,time) );
640
641 }
642
643 //________________________________________________________________________________
644 AliTRDarrayADC *AliTRDdigitsManager::GetDigits(Int_t det) const
645 {
646   //
647   // Returns the digits array for one detector
648   //
649
650   Int_t RecoDet = fRawRec ? 0 : det;
651
652   if (!fDigits)   
653     {
654       return 0x0;
655     }
656
657   if (!fHasSDigits)
658     {
659       ((AliTRDarrayADC *) fDigits->At(RecoDet))->SetNdet(det);
660       return (AliTRDarrayADC *) fDigits->At(RecoDet); 
661     }
662   else
663     {
664       AliDebug(2,"ERROR IN DATA TYPE!!!!");
665       return 0x0;
666     }
667
668 }
669
670 //_____________________________________________________________________________
671 AliTRDarraySignal *AliTRDdigitsManager::GetSDigits(Int_t det) const
672 {
673   //
674   // Returns the sdigits array for one detector
675   //
676
677   Int_t RecoDet = fRawRec ? 0 : det;
678
679   if (!fDigits)   
680     {
681       //      AliDebug(1,"NO FDIGITS!");        
682       return 0x0;
683     }
684
685   if (fHasSDigits)
686     {
687       ((AliTRDarraySignal *) fDigits->At(RecoDet))->SetNdet(det);
688       return (AliTRDarraySignal *) fDigits->At(RecoDet);
689     }
690   else
691     {
692       AliDebug(2,"ERROR IN DATA TYPE!!!!");
693       return 0x0;
694     }
695
696 }
697
698 //_____________________________________________________________________________
699 AliTRDarrayDictionary *AliTRDdigitsManager::GetDictionary(Int_t det
700                                                         , Int_t i) const
701 {
702   //
703   // Returns the dictionary for one detector
704   //
705
706   Int_t RecoDet = fRawRec ? 0 : det;
707
708   if (fUseDictionaries == kFALSE)
709     {
710       return 0x0;
711     }
712
713   ((AliTRDarrayDictionary *) fDigits->At(RecoDet))->SetNdet(det);
714   return (AliTRDarrayDictionary *) fDict[i]->At(RecoDet);
715   
716 }
717
718 //_____________________________________________________________________________
719 Int_t AliTRDdigitsManager::GetTrack(Int_t track, AliTRDdigit * const digit) const
720 {
721   // 
722   // Returns the MC-track numbers from the dictionary for a given digit
723   //
724
725   Int_t row  = digit->GetRow();
726   Int_t col  = digit->GetCol();
727   Int_t time = digit->GetTime();
728   Int_t det  = digit->GetDetector();
729
730   return GetTrack(track,row,col,time,det);
731
732 }
733
734 //_____________________________________________________________________________
735 AliTRDSignalIndex *AliTRDdigitsManager::GetIndexes(Int_t det) 
736 {
737   // 
738   // Returns indexes of active pads
739   //
740
741   Int_t RecoDet = fRawRec ? 0 : det;
742
743   return (AliTRDSignalIndex *) fSignalIndexes->At(RecoDet);
744
745 }
746
747 //_____________________________________________________________________________
748 void AliTRDdigitsManager::RemoveDigits(Int_t det) 
749 {
750    // 
751    // Clear memory at det for Digits
752    //
753
754   Int_t RecoDet = fRawRec ? 0 : det;
755
756   if (fDigits->At(RecoDet))
757     {
758       if (fHasSDigits) 
759         {
760           AliTRDarraySignal *arr = (AliTRDarraySignal *) fDigits->RemoveAt(RecoDet);
761           delete arr;
762         }
763       else 
764         {
765           AliTRDarrayADC    *arr = (AliTRDarrayADC *)    fDigits->RemoveAt(RecoDet);
766           delete arr;
767         }
768     }
769
770 }
771
772 //_____________________________________________________________________________
773 void AliTRDdigitsManager::RemoveDictionaries(Int_t det) 
774 {
775   // 
776   // Clear memory
777   //
778
779   Int_t RecoDet = fRawRec ? 0 : det;
780
781   if (fUseDictionaries == kFALSE) 
782     {
783       return;
784     }
785
786   for (Int_t i = 0; i < kNDict; i++) 
787     {
788       if (fDict[i]->At(RecoDet))
789         {
790           AliTRDarrayDictionary *arr = (AliTRDarrayDictionary *) fDict[i]->RemoveAt(RecoDet);
791           delete arr;
792         }
793     }
794
795 }
796
797 //_____________________________________________________________________________
798 void AliTRDdigitsManager::RemoveIndexes(Int_t det) 
799 {
800    // 
801    // Clear memory
802    //
803
804   Int_t RecoDet = fRawRec ? 0 : det;
805
806   if (fSignalIndexes->At(RecoDet))
807     {
808       AliTRDSignalIndex *arr = (AliTRDSignalIndex *) fSignalIndexes->RemoveAt(RecoDet);
809       delete arr;
810     }
811
812 }
813
814 //_____________________________________________________________________________
815 void AliTRDdigitsManager::ClearIndexes(Int_t det) 
816 {
817   // 
818   // Clear memory
819   //
820   
821   Int_t RecoDet = fRawRec ? 0 : det;
822
823   ((AliTRDSignalIndex *) fSignalIndexes->At(RecoDet))->ClearAll();  
824
825 }
826
827 //_____________________________________________________________________________
828 Bool_t AliTRDdigitsManager::BuildIndexes(Int_t det)
829 {
830   //
831   // Build the list of indices
832   //
833
834   Int_t nRows  = 0;
835   Int_t nCols  = 0;
836   Int_t nTbins = 0;
837
838   AliTRDgeometry  geom;
839   AliTRDarrayADC *digits = 0x0;
840
841   if (fHasSDigits) 
842     {
843       return kFALSE;
844     }
845   else 
846     {
847       digits = (AliTRDarrayADC *) GetDigits(det);
848     }
849
850   // digits should be expanded by now!!!
851   if (digits->GetNtime() > 0) 
852     {      
853       digits->Expand(); 
854       nRows  = digits->GetNrow();
855       nCols  = digits->GetNcol();
856       nTbins = digits->GetNtime();
857       
858       AliTRDSignalIndex *indexes = GetIndexes(det);
859       indexes->SetSM(geom.GetSector(det));
860       indexes->SetStack(geom.GetStack(det));
861       indexes->SetLayer(geom.GetLayer(det));
862       indexes->SetDetNumber(det);
863
864       if (indexes->IsAllocated() == kFALSE)
865         {
866           indexes->Allocate(nRows,nCols,nTbins);
867         }
868
869       for (Int_t ir = 0; ir < nRows; ir++) 
870         {
871           for (Int_t ic = 0; ic < nCols; ic++) 
872             {
873               for (Int_t it = 0; it < nTbins; it++)
874                 {         
875                   Int_t isig = digits->GetDataBits(ir,ic,it);
876                   if (isig > 0) 
877                     {
878                       indexes->AddIndexRC(ir,ic);           
879                     }
880                 } // tbins
881             } // cols
882         } // rows
883
884     } // if GetNtime
885   else 
886     {
887       return kFALSE;
888     }
889   
890   return kTRUE;
891
892 }
893
894 //_____________________________________________________________________________
895 Bool_t AliTRDdigitsManager::LoadArrayDigits()
896 {
897   //
898   // Loads the (s-)digits arrays for all detectors
899   //
900
901   if (!fTree) 
902     {
903       AliError("Digits tree is not defined\n");
904       return kFALSE;
905     }
906
907   Bool_t status = kTRUE;
908
909   // Get the branch
910   TBranch *branch = fTree->GetBranch("TRDdigits");
911   if (!branch) 
912     {
913       AliError("Branch TRDdigits is not defined\n");
914       return kFALSE;
915     }
916
917   // Loop through all detectors and read them from the tree
918   for (Int_t iDet = 0; iDet < fDets; iDet++) 
919     {
920       if (fHasSDigits)
921         {
922           AliTRDarraySignal *dataArray = (AliTRDarraySignal *) fDigits->At(iDet);
923           if (!dataArray) 
924             {
925               status = kFALSE;
926               break;    
927             }
928           branch->SetAddress(&dataArray);
929           branch->GetEntry(iDet);
930         }
931       else
932         {
933           AliTRDarrayADC    *dataArray = (AliTRDarrayADC *)    fDigits->At(iDet);
934           if (!dataArray) 
935             {
936               status = kFALSE;
937               break;
938             }
939           branch->SetAddress(&dataArray);
940           branch->GetEntry(iDet);
941         }
942     }
943
944   return status;
945
946 }
947
948 //________________________________________________________________________________________________
949 Bool_t AliTRDdigitsManager::LoadArrayDict()
950 {
951   //
952   // Loads dictionary arrays for all detectors
953   //
954
955   if (!fTree) 
956     {
957       AliError("Digits tree is not defined\n");
958       return kFALSE;
959     }
960
961   Bool_t status = kTRUE;
962
963   for (Int_t iDict = 0; iDict < kNDict; iDict++) 
964     {
965
966       // Get the branch
967       Char_t branchname[15];
968       sprintf(branchname,"TRDdictionary%d",iDict);
969       TBranch *branch = fTree->GetBranch(branchname);
970       if (!branch) 
971         {
972           AliError(Form("Branch %s is not defined\n",branchname));
973           return kFALSE;
974         }
975
976       // Loop through all detectors and read them from the tree
977       for (Int_t iDet = 0; iDet < fDets; iDet++) 
978         {
979           AliTRDarrayDictionary *dataArray = (AliTRDarrayDictionary *) fDict[iDict]->At(iDet);
980           if (!dataArray) 
981             {
982               status = kFALSE;
983               break;    
984             }
985           branch->SetAddress(&dataArray);
986           branch->GetEntry(iDet);
987         }
988
989     }
990
991   return status;
992
993 }
994
995 //_____________________________________________________________________________
996 Bool_t AliTRDdigitsManager::LoadDigitsParam()
997 {
998   //
999   // Loads the digits parameter object from the digits tree
1000   //
1001
1002   if (!fTree) 
1003     {
1004       AliError("Digits tree is not defined\n");
1005       return kFALSE;
1006     }
1007
1008   // Get the branch
1009   TBranch *branch = fTree->GetBranch("TRDdigitsParam");
1010   if (!branch) 
1011     {
1012       AliError("Branch TRDdigitsParam is not defined\n");
1013       return kFALSE;
1014     }
1015
1016   // Read the parameter object
1017   AliTRDdigitsParam *digitsParam = fDigitsParam;
1018   if (!digitsParam)
1019     {
1020       return kFALSE;
1021     }
1022
1023   branch->SetAddress(&digitsParam);
1024   branch->GetEntry();
1025
1026   return kTRUE;
1027
1028 }
1029
1030 //_____________________________________________________________________________
1031 Bool_t AliTRDdigitsManager::StoreArrayDigits()
1032 {
1033   //
1034   // Stores the digit arrays for all detectors
1035   //
1036
1037   if (!fTree) 
1038     {
1039       AliError("Digits tree is not defined\n");
1040       return kFALSE;
1041     }
1042
1043   // Get the branch
1044   TBranch *branch = fTree->GetBranch("TRDdigits");
1045   if (!branch) 
1046     {
1047       AliError("Branch TRDdigits is not defined\n");
1048       return kFALSE;
1049     }
1050
1051   // Loop through all detectors and fill them into the tree
1052   Bool_t status = kTRUE;
1053   for (Int_t iDet = 0; iDet < fDets; iDet++) 
1054     {
1055       if (fHasSDigits)
1056         {
1057           const AliTRDarraySignal *kDataArray = (AliTRDarraySignal *) fDigits->At(iDet);
1058           if (!kDataArray) 
1059             {
1060               status = kFALSE;
1061               break;
1062             }
1063           branch->SetAddress(&kDataArray);
1064           branch->Fill();
1065         }
1066       else
1067         {
1068           const AliTRDarrayADC    *kDataArray = (AliTRDarrayADC *)    fDigits->At(iDet); 
1069           if (!kDataArray) 
1070             {
1071               status = kFALSE;
1072               break;
1073             }
1074           branch->SetAddress(&kDataArray);
1075           branch->Fill();
1076         }
1077     }
1078
1079   return status;
1080
1081 }
1082
1083 //_____________________________________________________________________________
1084 Bool_t AliTRDdigitsManager::StoreArrayDict()
1085 {
1086   //
1087   // Stores the dictionary arrays for all detectors
1088   //
1089
1090   if (!fTree) 
1091     {
1092       AliError("Digits tree is not defined\n");
1093       return kFALSE;
1094     }
1095
1096   Bool_t status = kTRUE;
1097
1098   for (Int_t iDict = 0; iDict < kNDict; iDict++)
1099      {
1100
1101        // Get the branch
1102        Char_t branchname[15];
1103        sprintf(branchname,"TRDdictionary%d",iDict);
1104        TBranch *branch = fTree->GetBranch(branchname);
1105        if (!branch) 
1106          {
1107            AliError(Form("Branch %s is not defined\n",branchname));
1108            return kFALSE;
1109          }
1110
1111        // Loop through all detectors and fill them into the tree
1112        for (Int_t iDet = 0; iDet < fDets; iDet++) 
1113          {
1114            const AliTRDarrayDictionary *kDataArray = (AliTRDarrayDictionary *) fDict[iDict]->At(iDet);
1115            if (!kDataArray) 
1116              {
1117                status = kFALSE;
1118                break;
1119              }
1120            branch->SetAddress(&kDataArray);
1121            branch->Fill();
1122          }
1123
1124      }
1125
1126   return status;
1127
1128 }
1129
1130 //_____________________________________________________________________________
1131 Bool_t AliTRDdigitsManager::StoreDigitsParam()
1132 {
1133   //
1134   // Stores the digits parameter object from the digits tree
1135   //
1136
1137   if (!fTree) 
1138     {
1139       AliError("Digits tree is not defined\n");
1140       return kFALSE;
1141     }
1142
1143   // Get the branch
1144   TBranch *branch = fTree->GetBranch("TRDdigitsParam");
1145   if (!branch) 
1146     {
1147       AliError("Branch TRDdigitsParam is not defined\n");
1148       return kFALSE;
1149     }
1150
1151   // Fill the digits object in the tree
1152   const AliTRDdigitsParam *kDigitsParam = fDigitsParam;
1153   if (!kDigitsParam)
1154     {
1155       return kFALSE;
1156     }
1157   branch->SetAddress(&kDigitsParam);
1158   branch->Fill();
1159
1160   return kTRUE;
1161
1162 }
1163