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