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