]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRDsegmentArray.cxx
Memory leaks
[u/mrichter/AliRoot.git] / TRD / AliTRDsegmentArray.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 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 "AliLog.h"
32
33 #include "AliTRDgeometry.h"
34 #include "AliTRDsegmentArray.h"
35 #include "AliTRDsegmentID.h"
36 #include "AliTRDdataArray.h"
37 #include "AliTRDarrayI.h"
38
39 ClassImp(AliTRDsegmentArray)
40
41 //_____________________________________________________________________________
42 AliTRDsegmentArray::AliTRDsegmentArray()
43   :TNamed()
44   ,fSegment(0) 
45   ,fTreeIndex(0)
46   ,fNSegment(0)
47   ,fTree(0)
48   ,fBranch(0)
49   ,fClass(0)
50 {
51   //
52   // Default constructor
53   //
54
55 }
56
57 //_____________________________________________________________________________
58 AliTRDsegmentArray::AliTRDsegmentArray(const char *classname, Int_t n)
59   :TNamed()
60   ,fSegment(0) 
61   ,fTreeIndex(0)
62   ,fNSegment(0)
63   ,fTree(0)
64   ,fBranch(0)
65   ,fClass(0)
66 {
67   //
68   //  Create an array of objects of <classname>. The class must inherit from
69   //  AliTRDsegmentID. The second argument sets the number of entries in 
70   //  the array.
71   //
72
73   AliTRDdataArray *dataArray;  
74
75   SetClass(classname);
76
77   if (MakeArray(n) == kFALSE) {
78     AliError(Form("Cannot allocate %d segments in memory",n));
79     return;
80   }
81
82   for (Int_t i = 0; i < n; i++) {
83     dataArray = (AliTRDdataArray *) AddSegment(i);
84   }
85
86 }
87
88 //_____________________________________________________________________________
89 AliTRDsegmentArray::AliTRDsegmentArray(AliTRDsegmentArray &a)
90   :TNamed(a)
91   ,fSegment(a.fSegment) 
92   ,fTreeIndex(a.fTreeIndex)
93   ,fNSegment(a.fNSegment)
94   ,fTree(a.fTree)
95   ,fBranch(a.fBranch)
96   ,fClass(a.fClass)
97 {
98   //
99   // AliTRDsegmentArray copy constructor
100   //
101
102   a.Copy(*this);
103
104 }
105
106 //_____________________________________________________________________________
107 AliTRDsegmentArray::~AliTRDsegmentArray()
108 {
109   //
110   // AliTRDsegmentArray destructor
111   //
112
113   // Needed ????
114   //Delete();
115
116   if (fNSegment) {
117     fSegment->Delete();
118     delete fSegment;
119   }
120
121   if (fTreeIndex) {
122     delete fTreeIndex;
123   }
124
125 }
126
127 //_____________________________________________________________________________
128 AliTRDsegmentArray &AliTRDsegmentArray::operator=(const AliTRDsegmentArray &a)
129 {
130   //
131   // Assignment operator
132   //
133
134   if (this != &a) ((AliTRDsegmentArray &) a).Copy(*this);
135   return *this;
136
137 }
138
139 //_____________________________________________________________________________
140 void AliTRDsegmentArray::Copy(TObject &a) const
141 {
142   //
143   // Copy function
144   //
145
146   TNamed::Copy(a);
147
148   fSegment->Copy(*((AliTRDsegmentArray &) a).fSegment);
149   fTreeIndex->Copy(*((AliTRDsegmentArray &) a).fTreeIndex);
150   fClass->Copy(*((AliTRDsegmentArray &) a).fClass);
151
152   ((AliTRDsegmentArray &) a).fNSegment = fNSegment;
153
154 }
155
156 //_____________________________________________________________________________
157 Bool_t AliTRDsegmentArray::SetClass(const Char_t *classname)
158 {
159   //
160   // Sets the classname of the stored object
161   //
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) {
177     AliFatal("ROOT system not initialized");
178     exit(1);
179   }   
180
181   fClass = gROOT->GetClass(classname);
182   if (!fClass) {
183     AliError(Form("%s is not a valid class name",classname));
184     return kFALSE;
185   }
186   if (!fClass->InheritsFrom(AliTRDsegmentID::Class())) {
187     AliError(Form("%s does not inherit from AliTRDsegmentID",classname));
188     return kFALSE;
189   }
190   
191   return kTRUE;
192
193 }
194
195 //_____________________________________________________________________________
196 AliTRDsegmentID *AliTRDsegmentArray::NewSegment()
197 {
198   //
199   // Create a new object according to the class information
200   //
201
202   if (fClass  == 0) {
203     return 0;
204   }
205
206   AliTRDsegmentID *segment = (AliTRDsegmentID *) fClass->New();
207
208   if (segment == 0) {
209     return 0;
210   }
211   else {
212     return segment;
213   }
214
215 }
216
217 //_____________________________________________________________________________
218 Bool_t AliTRDsegmentArray::AddSegment(AliTRDsegmentID *segment)
219 {
220   //
221   // Add a segment to the array
222   //
223
224   if (segment  == 0) {
225     return kFALSE;
226   }
227   if (fSegment == 0) {
228     return kFALSE;
229   }
230   if (fClass   == 0) {
231     return kFALSE;
232   }
233
234   if (!(segment->IsA()->InheritsFrom(fClass))) {
235     AliError(Form("added class %s is not of proper type"
236                  ,segment->IsA()->GetName()));
237     return kFALSE;
238   }
239
240   fSegment->AddAt(segment,segment->GetID());
241   fNSegment = fSegment->GetLast() + 1;
242
243   return kTRUE;
244
245 }
246
247 //_____________________________________________________________________________
248 AliTRDsegmentID *AliTRDsegmentArray::AddSegment(Int_t index)
249 {
250   //
251   // Add a segment to the array
252   //
253
254   if (fSegment == 0) {
255     return 0;
256   }
257   if (fClass   == 0) {
258     return 0;
259   }
260
261   AliTRDsegmentID *segment = NewSegment();
262   if (segment  == 0) {
263     return 0;
264   }
265
266   fSegment->AddAt(segment,index);
267   segment->SetID(index);
268   fNSegment = fSegment->GetLast() + 1;
269
270   return segment;
271
272 }
273
274 //_____________________________________________________________________________
275 Bool_t AliTRDsegmentArray::MakeArray(Int_t n)
276 {
277   //
278   // Create an array of pointers to the segments
279   //
280
281   if (fSegment) {
282     fSegment->Delete();
283     delete fSegment;
284   }
285   if (fTreeIndex) delete fTreeIndex;  
286
287   fSegment   = new TObjArray(n);
288   fTreeIndex = new AliTRDarrayI();
289   fTreeIndex->Set(n);
290   fNSegment  = n;
291   if ((fSegment) && (fTreeIndex)) {
292     return kTRUE;
293   }
294   else { 
295     return kFALSE;
296   }
297                   
298 }
299
300 //_____________________________________________________________________________
301 void AliTRDsegmentArray::ClearSegment(Int_t index)
302 {
303   //
304   // Remove a segment from the active memory    
305   //
306
307   if (fSegment->At(index)) {
308     fClass->Destructor(fSegment->RemoveAt(index)); 
309   }
310
311 }
312
313 //_____________________________________________________________________________
314 void AliTRDsegmentArray::MakeTree(char *file)
315 {
316   //
317   // Create a tree for the segment
318   //
319
320   AliTRDsegmentID *psegment = NewSegment();  
321
322   if (fTree) {
323     delete fTree;
324   }
325
326   fTree   = new TTree("Segment Tree","Tree with segments");
327   fBranch = fTree->Branch("Segment",psegment->IsA()->GetName(),&psegment,64000);
328
329   if (file) {
330     fBranch->SetFile(file);      
331   }
332
333   delete psegment;
334
335 }              
336
337 //_____________________________________________________________________________
338 Bool_t AliTRDsegmentArray::ConnectTree(const char *treeName)
339 {
340   //
341   // Connect a tree from current directory  
342   //
343
344   if (fTree) {
345     delete fTree;
346     fTree   = 0;
347     fBranch = 0;
348   }
349
350   fTree = (TTree *) gDirectory->Get(treeName);
351   if (fTree   == 0) {
352     return kFALSE;
353   }
354   fBranch = fTree->GetBranch("Segment");
355   if (fBranch == 0) {
356     return kFALSE;
357   }
358
359   MakeDictionary(TMath::Max(fNSegment,Int_t(fTree->GetEntries())));
360
361   return kTRUE;
362
363 }
364
365 //_____________________________________________________________________________
366 AliTRDsegmentID *AliTRDsegmentArray::LoadSegment(Int_t index)
367 {
368   //
369   // Load a segment with index <index> into the memory
370   //
371
372   if (fTreeIndex == 0) {
373     MakeDictionary(3000);
374   }
375
376   // First try to load dictionary 
377   if (fTreeIndex == 0) {
378     return 0;
379   }
380   if (fBranch    == 0) {
381     return 0;
382   }
383   if (index > fTreeIndex->fN) {
384     return 0;
385   }
386
387   AliTRDsegmentID *s = (AliTRDsegmentID *) fSegment->At(index);
388   if (s == 0) {
389     s = NewSegment();
390   }
391   s->SetID(index);
392   
393   if (s != 0) {
394     Int_t treeIndex = (*fTreeIndex)[index];
395     if (treeIndex < 1) {
396       return 0;
397     }
398     else { 
399       treeIndex--;
400     }   
401     fBranch->SetAddress(&s);
402     fTree->GetEvent(treeIndex);
403     fSegment->AddAt((TObject*) s, index);
404   }
405   else { 
406     return 0;
407   }
408
409   return s;
410
411 }
412
413 //_____________________________________________________________________________
414 AliTRDsegmentID *AliTRDsegmentArray::LoadEntry(Int_t index)
415 {
416   //
417   // Load a segment at position <index> in the tree into the memory
418   //
419
420   if (fBranch == 0) {
421     return 0;
422   }
423   if (index > fTree->GetEntries()) {
424     return 0;
425   }
426
427   AliTRDsegmentID *s = NewSegment();  
428   if (s) {
429     fBranch->SetAddress(&s);
430     fTree->GetEvent(index);
431   }
432   else {
433     return 0;
434   }
435
436   Int_t nindex = s->GetID();
437   ClearSegment(nindex);
438   fSegment->AddAt((TObject *) s, nindex);
439
440   return s;
441
442 }
443
444 //_____________________________________________________________________________
445 void AliTRDsegmentArray::StoreSegment(Int_t index)
446 {
447   //
448   // Make a segment persistent 
449   //
450
451   const AliTRDsegmentID *kSegment = (*this)[index];
452   if (kSegment == 0) {
453     return;
454   }
455   if (fTree    == 0) {
456     MakeTree();
457   }
458   fBranch->SetAddress(&kSegment);
459   fTree->Fill();
460
461 }
462
463 //_____________________________________________________________________________
464 Bool_t  AliTRDsegmentArray::MakeDictionary(Int_t size)
465 {
466   //
467   // Create an index table for the tree
468   //  
469
470   if (size < 1) {
471     return kFALSE;
472   }
473   if (fTreeIndex) {
474     delete fTreeIndex;
475   }
476
477   fTreeIndex = new AliTRDarrayI(); 
478   fTreeIndex->Set(size);
479   
480   AliTRDsegmentID   segment;
481   AliTRDsegmentID *psegment = &segment;
482
483   fBranch->SetAddress(&psegment);
484   TBranch *brindix = fTree->GetBranch("fSegmentID");
485
486   Int_t nevent = (Int_t) fTree->GetEntries();  
487   for (Int_t i = 0; i < nevent; i++) {
488     brindix->GetEvent(i);
489     Int_t treeIndex = segment.GetID();
490     if (fTreeIndex->fN < treeIndex) {
491       fTreeIndex->Expand(Int_t (Float_t(treeIndex) * 1.5) + 1);
492     }
493     (*fTreeIndex)[treeIndex] = i + 1; 
494   }
495
496   return kTRUE;
497
498 }
499
500 //_____________________________________________________________________________
501 const AliTRDsegmentID * AliTRDsegmentArray::operator[](Int_t i) const
502 {
503   //
504   // Returns a segment with the given index <i>
505   //
506
507   if ((i <          0) || 
508       (i >= fNSegment)) {
509     return 0; 
510   }
511
512   return (AliTRDsegmentID *) fSegment->At(i);
513
514 }
515
516 //_____________________________________________________________________________
517 const AliTRDsegmentID *AliTRDsegmentArray::At(Int_t i) const
518 {
519   //
520   // Returns a segment with the given index <i>
521   //
522
523   if ((i <          0) || 
524       (i >= fNSegment)) {
525     return 0; 
526   }
527
528   return (AliTRDsegmentID *) fSegment->At(i);
529
530 }
531
532 //_____________________________________________________________________________
533 void AliTRDsegmentArray::Delete()
534 {
535   //
536   // Deletes all detector segments from the array
537   //
538
539   for (Int_t iDet = 0; iDet < fNSegment; iDet++) {
540     ClearSegment(iDet);
541   }
542
543 }
544
545 //_____________________________________________________________________________
546 Bool_t AliTRDsegmentArray::LoadArray(const Char_t *branchname, TTree *tree)
547 {
548   //
549   // Loads all segments of the array from the branch <branchname> of
550   // the digits tree <tree>
551   //
552
553   fTree = tree;
554
555   if (!fTree) {
556     AliError("Digits tree is not defined\n");
557     return kFALSE;
558   }
559
560   // Get the branch
561   fBranch = fTree->GetBranch(branchname);
562   if (!fBranch) {
563     AliError(Form("Branch %s is not defined\n",branchname));
564     return kFALSE;
565   }
566
567   // Loop through all segments and read them from the tree
568   Bool_t status = kTRUE;
569   for (Int_t iSegment = 0; iSegment < fNSegment; iSegment++) {
570     AliTRDdataArray *dataArray = (AliTRDdataArray *) fSegment->At(iSegment);
571     if (!dataArray) {
572       status = kFALSE;
573       break;    
574     }
575     fBranch->SetAddress(&dataArray);
576     fBranch->GetEntry(iSegment);
577   }
578
579   return status;
580
581 }
582
583 //_____________________________________________________________________________
584 Bool_t AliTRDsegmentArray::StoreArray(const Char_t *branchname, TTree *tree)
585 {
586   //
587   // Stores all segments of the array in the branch <branchname> of 
588   // the digits tree <tree>
589   //
590
591   fTree = tree;
592
593   if (!fTree) {
594     AliError("Digits tree is not defined\n");
595     return kFALSE;
596   }
597
598   // Get the branch
599   fBranch = fTree->GetBranch(branchname);
600   if (!fBranch) {
601     AliError(Form("Branch %s is not defined\n",branchname));
602     return kFALSE;
603   }
604
605   // Loop through all segments and fill them into the tree
606   Bool_t status = kTRUE;
607   for (Int_t iSegment = 0; iSegment < fNSegment; iSegment++) {
608     const AliTRDdataArray *kDataArray = 
609          (AliTRDdataArray *) AliTRDsegmentArray::At(iSegment);
610     if (!kDataArray) {
611       status = kFALSE;
612       break;
613     }
614     fBranch->SetAddress(&kDataArray);
615     fBranch->Fill();
616   }
617
618   return status;
619
620 }
621
622 //_____________________________________________________________________________
623 AliTRDdataArray *AliTRDsegmentArray::GetDataArray(Int_t det) const
624 {
625   //
626   // Returns the data array for a given detector
627   //
628
629   return ((AliTRDdataArray *) AliTRDsegmentArray::At(det));
630
631 }
632
633 //_____________________________________________________________________________
634 AliTRDdataArray *AliTRDsegmentArray::GetDataArray(Int_t pla
635                                                 , Int_t cha
636                                                 , Int_t sec) const
637 {
638   //
639   // Returns the data array for a given detector
640   //
641
642   Int_t det = AliTRDgeometry::GetDetector(pla,cha,sec);
643   return GetDataArray(det);
644
645 }