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