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