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