]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDdigitsManager.cxx
1976eca103c3c13a802c7a4282e0f308d5c69f95
[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   Int_t recoDet = fRawRec ? 0 : det;
283
284   RemoveDigits(recoDet);
285   RemoveDictionaries(recoDet);
286   RemoveIndexes(recoDet);
287
288   if (fHasSDigits)
289     fDigits->AddAt(new AliTRDarraySignal(),recoDet);
290   else
291     fDigits->AddAt(new AliTRDarrayADC(),recoDet);
292
293   if (fUseDictionaries) 
294     {
295       for (Int_t iDict = 0; iDict < kNDict; iDict++)
296         fDict[iDict]->AddAt(new AliTRDarrayDictionary(),recoDet);
297     }
298   
299   fSignalIndexes->AddAt(new AliTRDSignalIndex(),recoDet);
300 }
301
302 //_____________________________________________________________________________
303 Short_t AliTRDdigitsManager::GetDigitAmp(Int_t row, Int_t col,Int_t time, Int_t det) const
304 {
305   //
306   // Returns the amplitude of a digit
307   //
308
309   if (!GetDigits(det)) return 0;
310   
311   return ((Short_t) ((AliTRDarrayADC *) GetDigits(det))->GetDataBits(row,col,time));
312
313 }
314
315 //_____________________________________________________________________________
316 UChar_t AliTRDdigitsManager::GetPadStatus(Int_t row, Int_t col, Int_t time, Int_t det) const
317 {
318   //
319   // Returns the pad status for the requested pad
320   //
321         
322   if (!GetDigits(det)) return 0;
323
324   return ((UChar_t) ((AliTRDarrayADC *) GetDigits(det))->GetPadStatus(row,col,time));
325  
326 }
327
328 //_____________________________________________________________________________
329 Bool_t AliTRDdigitsManager::MakeBranch(TTree * const tree)  
330 {
331   //
332   // Creates the tree and branches for the digits and the dictionary
333   //
334
335   Int_t  buffersize = 64000;
336   Bool_t status     = kTRUE;
337
338   if (tree) 
339     {
340       fTree = tree;
341     }
342
343   // Make the branch for the digits
344   if (fDigits) 
345     {
346       if(fHasSDigits)
347         {
348           const AliTRDarraySignal *kDigits = (AliTRDarraySignal *) fDigits->At(0); 
349           if (kDigits) 
350             {
351               if (!fTree) return kFALSE;
352               AliDebug(1,"Making branch for SDigits!\n");
353               TBranch* branch = fTree->GetBranch("TRDdigits");
354               if (!branch) fTree->Branch("TRDdigits","AliTRDarraySignal",&kDigits,buffersize,99);
355               AliDebug(1,"Making branch TRDdigits\n");
356             }
357           else 
358             {
359               status = kFALSE;
360             }
361         }
362
363       if(!fHasSDigits)
364         {
365           const AliTRDarrayADC *kDigits = (AliTRDarrayADC *) fDigits->At(0);
366           if (kDigits) 
367             {
368               if (!fTree) return kFALSE;
369               AliDebug(1,"Making branch for Digits!\n");
370               TBranch* branch = fTree->GetBranch("TRDdigits");
371               if (!branch) fTree->Branch("TRDdigits","AliTRDarrayADC",&kDigits,buffersize,99);
372               AliDebug(1,"Making branch TRDdigits\n");        
373             }
374           else 
375             {
376               status = kFALSE;
377             }
378         }
379
380     }    
381   else
382     {
383       status = kFALSE;
384     }
385   
386   if (fUseDictionaries) 
387     {
388       // Make the branches for the dictionaries
389       for (Int_t iDict = 0; iDict < kNDict; iDict++) 
390         {
391           Char_t branchname[15];
392           sprintf(branchname,"TRDdictionary%d",iDict); 
393           if (fDict[iDict]) 
394             {
395               const AliTRDarrayDictionary *kDictionary = (AliTRDarrayDictionary *) fDict[iDict]->At(0);
396               if (kDictionary) 
397                 {
398                   if (!fTree) return kFALSE;
399                   AliDebug(2,"Making branch for dictionary!\n");
400                   TBranch* branch = fTree->GetBranch(branchname);
401                   if (!branch) fTree->Branch(branchname,"AliTRDarrayDictionary",&kDictionary,buffersize,99);
402                   AliDebug(1,Form("Making branch %s\n",branchname));
403                 }
404               else 
405                 {
406                   status = kFALSE;
407                 }
408             }
409           else 
410             {
411               status = kFALSE;
412             }
413         }
414     }
415   
416   return status;
417
418 }
419
420 //_____________________________________________________________________________
421 Bool_t AliTRDdigitsManager::ReadDigits(TTree * const tree)
422 {
423   //
424   // Reads the digit information from the input file
425   //
426
427   Bool_t status = kTRUE;
428
429   if (tree) 
430     {
431       fTree = tree;
432     }
433
434   if (!fDigits) 
435     {
436       AliDebug(1,"Create the data arrays.\n");
437       CreateArrays();
438     }
439
440   status = LoadArray(fDigits,"TRDdigits",fTree);
441
442   if (fUseDictionaries) 
443     {
444       for (Int_t iDict = 0; iDict < kNDict; iDict++) 
445         {
446           Char_t branchname[15];
447           sprintf(branchname,"TRDdictionary%d",iDict);
448           status = LoadArrayDict(fDict[iDict],branchname,fTree);
449           if (status == kFALSE) 
450             {
451               fUseDictionaries = kFALSE;
452               AliWarning("Unable to load dict arrays. Will not use them.\n");
453               break;
454             }
455         }  
456     }
457
458   return kTRUE;
459
460 }
461
462 //_____________________________________________________________________________
463 Bool_t AliTRDdigitsManager::WriteDigits()
464 {
465   //
466   // Writes out the TRD-digits and the dictionaries
467   //
468
469   // Store the contents of the detector array in the tree
470
471   if (!StoreArray(fDigits,"TRDdigits",fTree))
472     {
473       AliError("Error while storing digits in branch TRDdigits\n");
474       return kFALSE;
475     }
476
477   if (fUseDictionaries) 
478     {
479       for (Int_t iDict = 0; iDict < kNDict; iDict++)
480         {
481           Char_t branchname[15];
482           sprintf(branchname,"TRDdictionary%d",iDict);
483           if (!StoreArrayDict(fDict[iDict],branchname,fTree)) 
484             {
485               AliError(Form("Error while storing dictionary in branch %s\n",branchname));
486               return kFALSE;
487             }
488         }
489     }
490   
491   // Write the new tree to the output file
492   fTree->AutoSave();
493
494   return kTRUE;
495
496 }
497
498 //_____________________________________________________________________________
499 AliTRDdigit *AliTRDdigitsManager::GetDigit(Int_t row
500                                          , Int_t col
501                                          , Int_t time
502                                          , Int_t det) const
503 {
504   // 
505   // Creates a single digit object 
506   //
507
508   Int_t digits[4]; 
509   Int_t amp[1];
510
511   digits[0] = det;
512   digits[1] = row;
513   digits[2] = col;
514   digits[3] = time;
515
516   amp[0]    = ((AliTRDarrayADC *) GetDigits(det))->GetData(row,col,time);
517   
518   return (new AliTRDdigit(digits,amp));
519
520 }
521
522 //_____________________________________________________________________________
523 Int_t AliTRDdigitsManager::GetTrack(Int_t track
524                                   , Int_t row
525                                   , Int_t col
526                                   , Int_t time
527                                   , Int_t det) const
528 {
529   // 
530   // Returns the MC-track numbers from the dictionary.
531   //
532
533   if ((track < 0) || (track >= kNDict)) 
534     {
535       AliError(Form("track %d out of bounds (size: %d, this: 0x%08x)",track,kNDict,this));
536       return -1;
537     }
538
539   if (fUseDictionaries == kFALSE) 
540     {
541       return -1;
542     }
543
544   // Array contains index+1 to allow data compression--->Changed
545   return (((AliTRDarrayDictionary *) GetDictionary(det,track))->GetData(row,col,time) );
546
547 }
548
549 //________________________________________________________________________________
550 AliTRDarrayADC *AliTRDdigitsManager::GetDigits(Int_t det) const
551 {
552   //
553   // Returns the digits array for one detector
554   //
555
556   Int_t RecoDet = fRawRec ? 0 : det;
557
558   if (!fDigits)   
559     {
560       return 0x0;
561     }
562
563   if (!fHasSDigits)
564     {
565       ((AliTRDarrayADC *) fDigits->At(RecoDet))->SetNdet(det);
566       return (AliTRDarrayADC *) fDigits->At(RecoDet); 
567     }
568   else
569     {
570       AliDebug(2,"ERROR IN DATA TYPE!!!!");
571       return 0x0;
572     }
573
574 }
575
576 //_____________________________________________________________________________
577 AliTRDarraySignal *AliTRDdigitsManager::GetSDigits(Int_t det) const
578 {
579   //
580   // Returns the sdigits array for one detector
581   //
582
583   Int_t RecoDet = fRawRec ? 0 : det;
584
585   if (!fDigits)   
586     {
587       //      AliDebug(1,"NO FDIGITS!");        
588       return 0x0;
589     }
590
591   if (fHasSDigits)
592     {
593       ((AliTRDarraySignal *) fDigits->At(RecoDet))->SetNdet(det);
594       return (AliTRDarraySignal *) fDigits->At(RecoDet);
595     }
596   else
597     {
598       AliDebug(2,"ERROR IN DATA TYPE!!!!");
599       return 0x0;
600     }
601
602 }
603
604 //_____________________________________________________________________________
605 AliTRDarrayDictionary *AliTRDdigitsManager::GetDictionary(Int_t det
606                                                         , Int_t i) const
607 {
608   //
609   // Returns the dictionary for one detector
610   //
611
612   Int_t RecoDet = fRawRec ? 0 : det;
613
614   if (fUseDictionaries == kFALSE)
615     {
616       return 0x0;
617     }
618
619   ((AliTRDarrayDictionary *) fDigits->At(RecoDet))->SetNdet(det);
620   return (AliTRDarrayDictionary *) fDict[i]->At(RecoDet);
621   
622 }
623
624 //_____________________________________________________________________________
625 Int_t AliTRDdigitsManager::GetTrack(Int_t track, AliTRDdigit * const digit) const
626 {
627   // 
628   // Returns the MC-track numbers from the dictionary for a given digit
629   //
630
631   Int_t row  = digit->GetRow();
632   Int_t col  = digit->GetCol();
633   Int_t time = digit->GetTime();
634   Int_t det  = digit->GetDetector();
635
636   return GetTrack(track,row,col,time,det);
637
638 }
639
640 //_____________________________________________________________________________
641 AliTRDSignalIndex *AliTRDdigitsManager::GetIndexes(Int_t det) 
642 {
643   // 
644   // Returns indexes of active pads
645   //
646
647   Int_t RecoDet = fRawRec ? 0 : det;
648
649   return (AliTRDSignalIndex *) fSignalIndexes->At(RecoDet);
650
651 }
652
653 //_____________________________________________________________________________
654 void AliTRDdigitsManager::RemoveDigits(Int_t det) 
655 {
656    // 
657    // Clear memory at det for Digits
658    //
659
660   Int_t RecoDet = fRawRec ? 0 : det;
661
662   if (fDigits->At(RecoDet))
663     {
664       if (fHasSDigits) 
665         {
666           AliTRDarraySignal *arr = (AliTRDarraySignal *) fDigits->RemoveAt(RecoDet);
667           delete arr;
668         }
669       else 
670         {
671           AliTRDarrayADC    *arr = (AliTRDarrayADC *)    fDigits->RemoveAt(RecoDet);
672           delete arr;
673         }
674     }
675
676 }
677
678 //_____________________________________________________________________________
679 void AliTRDdigitsManager::RemoveDictionaries(Int_t det) 
680 {
681   // 
682   // Clear memory
683   //
684
685   Int_t RecoDet = fRawRec ? 0 : det;
686
687   if (fUseDictionaries == kFALSE) 
688     {
689       return;
690     }
691
692   for (Int_t i = 0; i < kNDict; i++) 
693     {
694       if (fDict[i]->At(RecoDet))
695         {
696           AliTRDarrayDictionary *arr = (AliTRDarrayDictionary *) fDict[i]->RemoveAt(RecoDet);
697           delete arr;
698         }
699     }
700
701 }
702
703 //_____________________________________________________________________________
704 void AliTRDdigitsManager::RemoveIndexes(Int_t det) 
705 {
706    // 
707    // Clear memory
708    //
709
710   Int_t RecoDet = fRawRec ? 0 : det;
711
712   if (fSignalIndexes->At(RecoDet))
713     {
714       AliTRDSignalIndex *arr = (AliTRDSignalIndex *) fSignalIndexes->RemoveAt(RecoDet);
715       delete arr;
716     }
717
718 }
719
720
721 //_____________________________________________________________________________
722 void AliTRDdigitsManager::ClearIndexes(Int_t det) 
723 {
724   // 
725   // Clear memory
726   //
727   
728   Int_t RecoDet = fRawRec ? 0 : det;
729
730   ((AliTRDSignalIndex *) fSignalIndexes->At(RecoDet))->ClearAll();  
731
732 }
733
734 //_____________________________________________________________________________
735 Bool_t AliTRDdigitsManager::BuildIndexes(Int_t det)
736 {
737   //
738   // Build the list of indices
739   //
740
741   Int_t nRows  = 0;
742   Int_t nCols  = 0;
743   Int_t nTbins = 0;
744
745   AliTRDgeometry  geom;
746   AliTRDarrayADC *digits = 0x0;
747
748   if (fHasSDigits) 
749     {
750       return kFALSE;
751     }
752   else 
753     {
754       digits = (AliTRDarrayADC *) GetDigits(det);
755     }
756
757   //digits should be expanded by now!!!
758   if (digits->GetNtime() > 0) 
759     {      
760       digits->Expand(); 
761       nRows  = digits->GetNrow();
762       nCols  = digits->GetNcol();
763       nTbins = digits->GetNtime();
764       
765       AliTRDSignalIndex *indexes = GetIndexes(det);
766       indexes->SetSM(geom.GetSector(det));
767       indexes->SetStack(geom.GetStack(det));
768       indexes->SetLayer(geom.GetLayer(det));
769       indexes->SetDetNumber(det);
770
771       if (indexes->IsAllocated() == kFALSE)
772         {
773           indexes->Allocate(nRows,nCols,nTbins);
774         }
775
776       for (Int_t ir = 0; ir < nRows; ir++) 
777         {
778           for (Int_t ic = 0; ic < nCols; ic++) 
779             {
780               for (Int_t it = 0; it < nTbins; it++)
781                 {         
782                   Int_t isig = digits->GetDataBits(ir,ic,it);
783                   if (isig > 0) 
784                     {
785                       indexes->AddIndexRC(ir,ic);           
786                     }
787                 } // tbins
788             } // cols
789         } // rows
790
791     } // if GetNtime
792   else 
793     {
794       return kFALSE;
795     }
796   
797   return kTRUE;
798
799 }
800
801 //_____________________________________________________________________________
802 Bool_t AliTRDdigitsManager::LoadArray(TObjArray * const object
803                                     , const Char_t *branchname
804                                     , TTree * const tree)
805 {
806   //
807   // Loads all detectors of the array from the branch <branchname> of
808   // the digits tree <tree>
809   // Adapted from code of the class AliTRDsegmentArray
810   //
811
812   fTreeD = tree;
813
814   if (!fTreeD) 
815     {
816       AliError("Digits tree is not defined\n");
817       return kFALSE;
818     }
819
820   // Get the branch
821   fBranch = fTreeD->GetBranch(branchname);
822   if (!fBranch) 
823     {
824       AliError(Form("Branch %s is not defined\n",branchname));
825       return kFALSE;
826     }
827
828   // Loop through all detectors and read them from the tree
829   Bool_t status = kTRUE;
830   for (Int_t iDet = 0; iDet < fDets; iDet++) 
831     {
832       if(fHasSDigits)
833         {
834           AliTRDarraySignal *dataArray = (AliTRDarraySignal *) object->At(iDet);
835           if (!dataArray) 
836             {
837               status = kFALSE;
838               break;    
839             }
840
841           fBranch->SetAddress(&dataArray);
842           fBranch->GetEntry(iDet);
843         }
844       else
845         {
846           AliTRDarrayADC *dataArray = (AliTRDarrayADC *) object->At(iDet);
847           if (!dataArray) 
848             {
849               status = kFALSE;
850               break;    
851             }
852           fBranch->SetAddress(&dataArray);
853           fBranch->GetEntry(iDet);
854         }
855     }
856
857   return status;
858
859 }
860
861 //________________________________________________________________________________________________
862 Bool_t AliTRDdigitsManager::LoadArrayDict(TObjArray * const object
863                                         , const Char_t *branchname
864                                         , TTree * const tree)
865 {
866   //
867   // Loads all detectors of the array from the branch <branchname> of
868   // the dictionary tree <tree>
869   // Adapted from code of the class AliTRDsegmentArray
870   //
871
872   fTreeD = tree;
873
874   if (!fTreeD) 
875     {
876       AliError("Digits tree is not defined\n");
877       return kFALSE;
878     }
879
880   // Get the branch
881   fBranch = fTreeD->GetBranch(branchname);
882   if (!fBranch) 
883     {
884       AliError(Form("Branch %s is not defined\n",branchname));
885       return kFALSE;
886     }
887
888   // Loop through all detectors and read them from the tree
889   Bool_t status = kTRUE;
890   for (Int_t iDet = 0; iDet < fDets; iDet++) 
891     {
892       AliTRDarrayDictionary *dataArray = (AliTRDarrayDictionary *) object->At(iDet);
893       if (!dataArray) 
894         {
895           status = kFALSE;
896           break;    
897         }
898       fBranch->SetAddress(&dataArray);
899       fBranch->GetEntry(iDet);
900     }
901
902   return status;
903
904 }
905
906 //_____________________________________________________________________________
907 Bool_t AliTRDdigitsManager::StoreArray(TObjArray * const array1
908                                      , const Char_t *branchname
909                                      , TTree * const tree)
910 {
911   //
912   // Stores all the detectors of the array in the branch <branchname> of 
913   // the digits tree <tree>
914   // Adapted from code of the class AliTRDsegmentArray
915   //
916
917   fTree = tree;
918
919   if (!fTree) 
920     {
921       AliError("Digits tree is not defined\n");
922       return kFALSE;
923     }
924
925   // Get the branch
926   fBranch = fTree->GetBranch(branchname);
927   if (!fBranch) 
928     {
929       AliError(Form("Branch %s is not defined\n",branchname));
930       return kFALSE;
931     }
932
933   // Loop through all detectors and fill them into the tree
934   Bool_t status = kTRUE;
935   for (Int_t iDet = 0; iDet < fDets; iDet++) 
936     {
937       if (fHasSDigits)
938         {
939           const AliTRDarraySignal *kDataArray = (AliTRDarraySignal *) array1->At(iDet);
940           if (!kDataArray) 
941             {
942               status = kFALSE;
943               break;
944             }
945           fBranch->SetAddress(&kDataArray);
946           fBranch->Fill();
947         }
948       else
949         {
950           const AliTRDarrayADC *kDataArray = (AliTRDarrayADC *) array1->At(iDet); 
951           if (!kDataArray) 
952             {
953               status = kFALSE;
954               break;
955             }
956           fBranch->SetAddress(&kDataArray);
957           fBranch->Fill();
958         }
959     }
960
961   return status;
962
963 }
964
965 //_____________________________________________________________________________
966 Bool_t AliTRDdigitsManager::StoreArrayDict(TObjArray * const array3
967                                          , const Char_t *branchname
968                                          , TTree * const tree)
969 {
970   //
971   // Stores all the dictionary arrays of the detectors of the array in the branch <branchname> of 
972   // the dictionary tree <tree>
973   // Adapted from code of the class AliTRDsegmentArray
974   //
975
976   //  AliDebug(1,"Storing Arrays of Dictionary");
977   fTree = tree;
978
979   if (!fTree) 
980     {
981       AliError("Digits tree is not defined\n");
982       return kFALSE;
983     }
984
985   // Get the branch
986   fBranch = fTree->GetBranch(branchname);
987   if (!fBranch) 
988     {
989       AliError(Form("Branch %s is not defined\n",branchname));
990       return kFALSE;
991     }
992
993   // Loop through all detectors and fill them into the tree
994   Bool_t status = kTRUE;
995   for (Int_t iDet = 0; iDet < fDets; iDet++) 
996     {
997       const AliTRDarrayDictionary *kDataArray = (AliTRDarrayDictionary *) array3->At(iDet);
998       if (!kDataArray) 
999         {
1000           status = kFALSE;
1001           break;
1002         }
1003       fBranch->SetAddress(&kDataArray);
1004       fBranch->Fill();
1005     }
1006
1007   return status;
1008
1009 }