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