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