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