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