]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDsegmentArrayBase.cxx
Re-structured digitizer
[u/mrichter/AliRoot.git] / TRD / AliTRDsegmentArrayBase.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 //  Alice segment manager base class                                         //
21 //                                                                           //
22 ///////////////////////////////////////////////////////////////////////////////
23
24 #include <TROOT.h>
25 #include <TTree.h>
26 #include <TClonesArray.h>
27 #include <TDirectory.h>
28 #include <TError.h>
29 #include <TClass.h>
30
31 #include "AliLog.h"
32
33 #include "AliTRDsegmentArrayBase.h"
34 #include "AliTRDarrayI.h"
35 #include "AliTRDsegmentID.h"
36
37 ClassImp(AliTRDsegmentArrayBase)
38   
39 //_____________________________________________________________________________
40 AliTRDsegmentArrayBase::AliTRDsegmentArrayBase()
41   :TNamed()
42   ,fSegment(0) 
43   ,fTreeIndex(0)
44   ,fNSegment(0)
45   ,fTree(0)
46   ,fBranch(0)
47   ,fClass(0)
48 {
49   //
50   // AliTRDsegmentArrayBase default constructor
51   //
52
53 }
54
55 //_____________________________________________________________________________
56 AliTRDsegmentArrayBase::AliTRDsegmentArrayBase(const char *classname, Int_t n)
57   :TNamed()
58   ,fSegment(0) 
59   ,fTreeIndex(0)
60   ,fNSegment(0)
61   ,fTree(0)
62   ,fBranch(0)
63   ,fClass(0)
64 {
65   //
66   //  Create an array of objects of <classname>. The class must inherit from
67   //  AliTRDsegmentID. The second argument sets the number of entries in 
68   //  the array.
69   //
70
71   SetClass(classname);
72
73   if (MakeArray(n) == kFALSE) {
74     Error("AliTRDsegmentArrayBase","Cannot allocate %d segments in memory",n);
75     return;
76   }
77
78 }
79
80 //_____________________________________________________________________________
81 AliTRDsegmentArrayBase::AliTRDsegmentArrayBase(const AliTRDsegmentArrayBase &a)
82   :TNamed(a)
83   ,fSegment(a.fSegment) 
84   ,fTreeIndex(a.fTreeIndex)
85   ,fNSegment(a.fNSegment)
86   ,fTree(a.fTree)
87   ,fBranch(a.fBranch)
88   ,fClass(a.fClass)
89 {
90   //
91   // AliTRDsegmentArrayBase copy constructor
92   //
93
94 }
95
96 //_____________________________________________________________________________
97 AliTRDsegmentArrayBase::~AliTRDsegmentArrayBase()
98 {
99   //
100   // AliTRDsegmentArrayBase destructor
101   //
102
103   if (fNSegment) {
104     fSegment->Delete();
105     delete fSegment;
106   }
107
108   if (fTreeIndex) {
109     delete fTreeIndex;
110   }
111
112 }
113
114 //_____________________________________________________________________________
115 AliTRDsegmentArrayBase &AliTRDsegmentArrayBase
116                         ::operator=(const AliTRDsegmentArrayBase &a)
117 {
118   //
119   // Assignment operator
120   //
121
122   if (this != &a) ((AliTRDsegmentArrayBase &) a).Copy(*this);
123   return *this;
124
125 }
126
127 //_____________________________________________________________________________
128 void AliTRDsegmentArrayBase::Copy(TObject &a) const
129 {
130   //
131   // Copy function
132   //
133
134   TNamed::Copy(a);
135
136   fSegment->Copy(*((AliTRDsegmentArrayBase &) a).fSegment);
137   fTreeIndex->Copy(*((AliTRDsegmentArrayBase &) a).fTreeIndex);
138   fClass->Copy(*((AliTRDsegmentArrayBase &) a).fClass);
139
140   ((AliTRDsegmentArrayBase &) a).fNSegment = fNSegment;
141
142 }
143
144 //_____________________________________________________________________________
145 Bool_t AliTRDsegmentArrayBase::SetClass(const char *classname)
146 {
147   //
148   // Sets the classname of the stored object
149   //
150
151   if (fTree    != 0) {
152     delete fTree;
153     fTree      = 0;
154     fBranch    = 0;
155     delete fTreeIndex;
156     fTreeIndex = 0;
157   } 
158   if (fSegment != 0) {
159     fSegment->Delete();
160     delete fSegment;
161     fSegment   = 0;
162   }
163
164   if (!gROOT) {
165     AliFatal("ROOT system not initialized");
166     exit(1);
167   }   
168
169   fClass = gROOT->GetClass(classname);
170   if (!fClass) {
171     AliError(Form("%s is not a valid class name",classname));
172     return kFALSE;
173   }
174   if (!fClass->InheritsFrom(AliTRDsegmentID::Class())) {
175     AliError(Form("%s does not inherit from AliTRDsegmentID",classname));
176     return kFALSE;
177   }
178   
179   return kTRUE;
180
181 }
182
183 //_____________________________________________________________________________
184 AliTRDsegmentID *AliTRDsegmentArrayBase::NewSegment()
185 {
186   //
187   // Create a new object according to the class information
188   //
189
190   if (fClass  == 0) {
191     return 0;
192   }
193
194   AliTRDsegmentID *segment = (AliTRDsegmentID *) fClass->New();
195
196   if (segment == 0) {
197     return 0;
198   }
199   else {
200     return segment;
201   }
202
203 }
204
205 //_____________________________________________________________________________
206 Bool_t AliTRDsegmentArrayBase::AddSegment(AliTRDsegmentID *segment)
207 {
208   //
209   // Add a segment to the array
210   //
211
212   if (segment  == 0) {
213     return kFALSE;
214   }
215   if (fSegment == 0) {
216     return kFALSE;
217   }
218   if (fClass   == 0) {
219     return kFALSE;
220   }
221
222   if (!(segment->IsA()->InheritsFrom(fClass))) {
223     AliError(Form("added class %s is not of proper type"
224                  ,segment->IsA()->GetName()));
225     return kFALSE;
226   }
227
228   fSegment->AddAt(segment,segment->GetID());
229   fNSegment = fSegment->GetLast() + 1;
230
231   return kTRUE;
232
233 }
234
235 //_____________________________________________________________________________
236 AliTRDsegmentID *AliTRDsegmentArrayBase::AddSegment(Int_t index)
237 {
238   //
239   // Add a segment to the array
240   //
241
242   if (fSegment == 0) {
243     return 0;
244   }
245   if (fClass   == 0) {
246     return 0;
247   }
248
249   AliTRDsegmentID *segment = NewSegment();
250   if (segment  == 0) {
251     return 0;
252   }
253
254   fSegment->AddAt(segment,index);
255   segment->SetID(index);
256   fNSegment = fSegment->GetLast() + 1;
257
258   return segment;
259
260 }
261
262 //_____________________________________________________________________________
263 Bool_t AliTRDsegmentArrayBase::MakeArray(Int_t n)
264 {
265   //
266   // Create an array of pointers to the segments
267   //
268
269   if (fSegment) {
270     fSegment->Delete();
271     delete fSegment;
272   }
273   if (fTreeIndex) delete fTreeIndex;  
274
275   fSegment   = new TObjArray(n);
276   fTreeIndex = new AliTRDarrayI();
277   fTreeIndex->Set(n);
278   fNSegment  = n;
279   if ((fSegment) && (fTreeIndex)) {
280     return kTRUE;
281   }
282   else { 
283     return kFALSE;
284   }
285                   
286 }
287
288 //_____________________________________________________________________________
289 void AliTRDsegmentArrayBase::ClearSegment(Int_t index)
290 {
291   //
292   // Remove a segment from the active memory    
293   //
294
295   if (fSegment->At(index)) {
296     delete fSegment->RemoveAt(index);
297   }
298
299 }
300
301 //_____________________________________________________________________________
302 void AliTRDsegmentArrayBase::MakeTree(char *file)
303 {
304   //
305   // Create a tree for the segment
306   //
307
308   AliTRDsegmentID *psegment = NewSegment();  
309
310   if (fTree) {
311     delete fTree;
312   }
313
314   fTree   = new TTree("Segment Tree","Tree with segments");
315   fBranch = fTree->Branch("Segment",psegment->IsA()->GetName(),&psegment,64000);
316
317   if (file) {
318     fBranch->SetFile(file);      
319   }
320
321   delete psegment;
322
323 }              
324
325 //_____________________________________________________________________________
326 Bool_t AliTRDsegmentArrayBase::ConnectTree(const char *treeName)
327 {
328   //
329   // Connect a tree from current directory  
330   //
331
332   if (fTree) {
333     delete fTree;
334     fTree   = 0;
335     fBranch = 0;
336   }
337
338   fTree = (TTree *) gDirectory->Get(treeName);
339   if (fTree   == 0) {
340     return kFALSE;
341   }
342   fBranch = fTree->GetBranch("Segment");
343   if (fBranch == 0) {
344     return kFALSE;
345   }
346
347   MakeDictionary(TMath::Max(fNSegment,Int_t(fTree->GetEntries())));
348
349   return kTRUE;
350
351 }
352
353 //_____________________________________________________________________________
354 AliTRDsegmentID *AliTRDsegmentArrayBase::LoadSegment(Int_t index)
355 {
356   //
357   // Load a segment with index <index> into the memory
358   //
359
360   if (fTreeIndex == 0) {
361     MakeDictionary(3000);
362   }
363
364   // First try to load dictionary 
365   if (fTreeIndex == 0) {
366     return 0;
367   }
368   if (fBranch    == 0) {
369     return 0;
370   }
371   if (index > fTreeIndex->fN) {
372     return 0;
373   }
374
375   AliTRDsegmentID *s = (AliTRDsegmentID *) fSegment->At(index);
376   if (s == 0) {
377     s = NewSegment();
378   }
379   s->SetID(index);
380   
381   if (s != 0) {
382     Int_t treeIndex = (*fTreeIndex)[index];
383     if (treeIndex < 1) {
384       return 0;
385     }
386     else { 
387       treeIndex--;
388     }   
389     fBranch->SetAddress(&s);
390     fTree->GetEvent(treeIndex);
391     fSegment->AddAt((TObject*) s, index);
392   }
393   else { 
394     return 0;
395   }
396
397   return s;
398
399 }
400
401 //_____________________________________________________________________________
402 AliTRDsegmentID *AliTRDsegmentArrayBase::LoadEntry(Int_t index)
403 {
404   //
405   // Load a segment at position <index> in the tree into the memory
406   //
407
408   if (fBranch == 0) {
409     return 0;
410   }
411   if (index > fTree->GetEntries()) {
412     return 0;
413   }
414
415   AliTRDsegmentID *s = NewSegment();  
416   if (s) {
417     fBranch->SetAddress(&s);
418     fTree->GetEvent(index);
419   }
420   else {
421     return 0;
422   }
423
424   Int_t nindex = s->GetID();
425   ClearSegment(nindex);
426   fSegment->AddAt((TObject *) s, nindex);
427
428   return s;
429
430 }
431
432 //_____________________________________________________________________________
433 void AliTRDsegmentArrayBase::StoreSegment(Int_t index)
434 {
435   //
436   // Make a segment persistent 
437   //
438
439   const AliTRDsegmentID *kSegment = (*this)[index];
440   if (kSegment == 0) {
441     return;
442   }
443   if (fTree    == 0) {
444     MakeTree();
445   }
446   fBranch->SetAddress(&kSegment);
447   fTree->Fill();
448
449 }
450
451 //_____________________________________________________________________________
452 Bool_t  AliTRDsegmentArrayBase::MakeDictionary(Int_t size)
453 {
454   //
455   // Create an index table for the tree
456   //  
457
458   if (size < 1) {
459     return kFALSE;
460   }
461   if (fTreeIndex) {
462     delete fTreeIndex;
463   }
464
465   fTreeIndex = new AliTRDarrayI(); 
466   fTreeIndex->Set(size);
467   
468   AliTRDsegmentID   segment;
469   AliTRDsegmentID *psegment = &segment;
470
471   fBranch->SetAddress(&psegment);
472   TBranch *brindix = fTree->GetBranch("fSegmentID");
473
474   Int_t nevent = (Int_t) fTree->GetEntries();  
475   for (Int_t i = 0; i < nevent; i++) {
476     brindix->GetEvent(i);
477     Int_t treeIndex = segment.GetID();
478     if (fTreeIndex->fN < treeIndex) {
479       fTreeIndex->Expand(Int_t (Float_t(treeIndex) * 1.5) + 1);
480     }
481     (*fTreeIndex)[treeIndex] = i + 1; 
482   }
483
484   return kTRUE;
485
486 }
487
488 //_____________________________________________________________________________
489 const AliTRDsegmentID * AliTRDsegmentArrayBase::operator[](Int_t i) const
490 {
491   //
492   // Returns a segment with the given index <i>
493   //
494
495   if ((i <          0) || 
496       (i >= fNSegment)) {
497     return 0; 
498   }
499
500   return (AliTRDsegmentID *) fSegment->At(i);
501
502 }
503
504 //_____________________________________________________________________________
505 const AliTRDsegmentID *AliTRDsegmentArrayBase::At(Int_t i) const
506 {
507   //
508   // Returns a segment with the given index <i>
509   //
510
511   if ((i <          0) || 
512       (i >= fNSegment)) {
513     return 0; 
514   }
515
516   return (AliTRDsegmentID *) fSegment->At(i);
517
518 }