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