1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 ////////////////////////////////////////////////////////////////////////////
20 // Alice segment manager class //
22 ////////////////////////////////////////////////////////////////////////////
26 #include <TClonesArray.h>
27 #include <TDirectory.h>
33 #include "AliTRDgeometry.h"
34 #include "AliTRDsegmentArray.h"
35 #include "AliTRDsegmentID.h"
36 #include "AliTRDdataArray.h"
37 #include "AliTRDarrayI.h"
39 ClassImp(AliTRDsegmentArray)
41 //_____________________________________________________________________________
42 AliTRDsegmentArray::AliTRDsegmentArray()
52 // Default constructor
57 //_____________________________________________________________________________
58 AliTRDsegmentArray::AliTRDsegmentArray(const char *classname, Int_t n)
68 // Create an array of objects of <classname>. The class must inherit from
69 // AliTRDsegmentID. The second argument sets the number of entries in
73 AliTRDdataArray *dataArray;
77 if (MakeArray(n) == kFALSE) {
78 AliError(Form("Cannot allocate %d segments in memory",n));
82 for (Int_t i = 0; i < n; i++) {
83 dataArray = (AliTRDdataArray *) AddSegment(i);
88 //_____________________________________________________________________________
89 AliTRDsegmentArray::AliTRDsegmentArray(AliTRDsegmentArray &a)
92 ,fTreeIndex(a.fTreeIndex)
93 ,fNSegment(a.fNSegment)
99 // AliTRDsegmentArray copy constructor
106 //_____________________________________________________________________________
107 AliTRDsegmentArray::~AliTRDsegmentArray()
110 // AliTRDsegmentArray destructor
127 //_____________________________________________________________________________
128 AliTRDsegmentArray &AliTRDsegmentArray::operator=(const AliTRDsegmentArray &a)
131 // Assignment operator
134 if (this != &a) ((AliTRDsegmentArray &) a).Copy(*this);
139 //_____________________________________________________________________________
140 void AliTRDsegmentArray::Copy(TObject &a) const
148 fSegment->Copy(*((AliTRDsegmentArray &) a).fSegment);
149 fTreeIndex->Copy(*((AliTRDsegmentArray &) a).fTreeIndex);
150 fClass->Copy(*((AliTRDsegmentArray &) a).fClass);
152 ((AliTRDsegmentArray &) a).fNSegment = fNSegment;
156 //_____________________________________________________________________________
157 Bool_t AliTRDsegmentArray::SetClass(const Char_t *classname)
160 // Sets the classname of the stored object
177 AliFatal("ROOT system not initialized");
181 fClass = gROOT->GetClass(classname);
183 AliError(Form("%s is not a valid class name",classname));
186 if (!fClass->InheritsFrom(AliTRDsegmentID::Class())) {
187 AliError(Form("%s does not inherit from AliTRDsegmentID",classname));
195 //_____________________________________________________________________________
196 AliTRDsegmentID *AliTRDsegmentArray::NewSegment()
199 // Create a new object according to the class information
206 AliTRDsegmentID *segment = (AliTRDsegmentID *) fClass->New();
217 //_____________________________________________________________________________
218 Bool_t AliTRDsegmentArray::AddSegment(AliTRDsegmentID *segment)
221 // Add a segment to the array
234 if (!(segment->IsA()->InheritsFrom(fClass))) {
235 AliError(Form("added class %s is not of proper type"
236 ,segment->IsA()->GetName()));
240 fSegment->AddAt(segment,segment->GetID());
241 fNSegment = fSegment->GetLast() + 1;
247 //_____________________________________________________________________________
248 AliTRDsegmentID *AliTRDsegmentArray::AddSegment(Int_t index)
251 // Add a segment to the array
261 AliTRDsegmentID *segment = NewSegment();
266 fSegment->AddAt(segment,index);
267 segment->SetID(index);
268 fNSegment = fSegment->GetLast() + 1;
274 //_____________________________________________________________________________
275 Bool_t AliTRDsegmentArray::MakeArray(Int_t n)
278 // Create an array of pointers to the segments
285 if (fTreeIndex) delete fTreeIndex;
287 fSegment = new TObjArray(n);
288 fTreeIndex = new AliTRDarrayI();
291 if ((fSegment) && (fTreeIndex)) {
300 //_____________________________________________________________________________
301 void AliTRDsegmentArray::ClearSegment(Int_t index)
304 // Remove a segment from the active memory
307 if (fSegment->At(index)) {
308 fClass->Destructor(fSegment->RemoveAt(index));
313 //_____________________________________________________________________________
314 void AliTRDsegmentArray::MakeTree(char *file)
317 // Create a tree for the segment
320 AliTRDsegmentID *psegment = NewSegment();
326 fTree = new TTree("Segment Tree","Tree with segments");
327 fBranch = fTree->Branch("Segment",psegment->IsA()->GetName(),&psegment,64000);
330 fBranch->SetFile(file);
337 //_____________________________________________________________________________
338 Bool_t AliTRDsegmentArray::ConnectTree(const char *treeName)
341 // Connect a tree from current directory
350 fTree = (TTree *) gDirectory->Get(treeName);
354 fBranch = fTree->GetBranch("Segment");
359 MakeDictionary(TMath::Max(fNSegment,Int_t(fTree->GetEntries())));
365 //_____________________________________________________________________________
366 AliTRDsegmentID *AliTRDsegmentArray::LoadSegment(Int_t index)
369 // Load a segment with index <index> into the memory
372 if (fTreeIndex == 0) {
373 MakeDictionary(3000);
376 // First try to load dictionary
377 if (fTreeIndex == 0) {
383 if (index > fTreeIndex->fN) {
387 AliTRDsegmentID *s = (AliTRDsegmentID *) fSegment->At(index);
394 Int_t treeIndex = (*fTreeIndex)[index];
401 fBranch->SetAddress(&s);
402 fTree->GetEvent(treeIndex);
403 fSegment->AddAt((TObject*) s, index);
413 //_____________________________________________________________________________
414 AliTRDsegmentID *AliTRDsegmentArray::LoadEntry(Int_t index)
417 // Load a segment at position <index> in the tree into the memory
423 if (index > fTree->GetEntries()) {
427 AliTRDsegmentID *s = NewSegment();
429 fBranch->SetAddress(&s);
430 fTree->GetEvent(index);
436 Int_t nindex = s->GetID();
437 ClearSegment(nindex);
438 fSegment->AddAt((TObject *) s, nindex);
444 //_____________________________________________________________________________
445 void AliTRDsegmentArray::StoreSegment(Int_t index)
448 // Make a segment persistent
451 const AliTRDsegmentID *kSegment = (*this)[index];
458 fBranch->SetAddress(&kSegment);
463 //_____________________________________________________________________________
464 Bool_t AliTRDsegmentArray::MakeDictionary(Int_t size)
467 // Create an index table for the tree
477 fTreeIndex = new AliTRDarrayI();
478 fTreeIndex->Set(size);
480 AliTRDsegmentID segment;
481 AliTRDsegmentID *psegment = &segment;
483 fBranch->SetAddress(&psegment);
484 TBranch *brindix = fTree->GetBranch("fSegmentID");
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);
493 (*fTreeIndex)[treeIndex] = i + 1;
500 //_____________________________________________________________________________
501 const AliTRDsegmentID * AliTRDsegmentArray::operator[](Int_t i) const
504 // Returns a segment with the given index <i>
512 return (AliTRDsegmentID *) fSegment->At(i);
516 //_____________________________________________________________________________
517 const AliTRDsegmentID *AliTRDsegmentArray::At(Int_t i) const
520 // Returns a segment with the given index <i>
528 return (AliTRDsegmentID *) fSegment->At(i);
532 //_____________________________________________________________________________
533 void AliTRDsegmentArray::Delete()
536 // Deletes all detector segments from the array
539 for (Int_t iDet = 0; iDet < fNSegment; iDet++) {
545 //_____________________________________________________________________________
546 Bool_t AliTRDsegmentArray::LoadArray(const Char_t *branchname, TTree *tree)
549 // Loads all segments of the array from the branch <branchname> of
550 // the digits tree <tree>
556 AliError("Digits tree is not defined\n");
561 fBranch = fTree->GetBranch(branchname);
563 AliError(Form("Branch %s is not defined\n",branchname));
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);
575 fBranch->SetAddress(&dataArray);
576 fBranch->GetEntry(iSegment);
583 //_____________________________________________________________________________
584 Bool_t AliTRDsegmentArray::StoreArray(const Char_t *branchname, TTree *tree)
587 // Stores all segments of the array in the branch <branchname> of
588 // the digits tree <tree>
594 AliError("Digits tree is not defined\n");
599 fBranch = fTree->GetBranch(branchname);
601 AliError(Form("Branch %s is not defined\n",branchname));
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);
614 fBranch->SetAddress(&kDataArray);
622 //_____________________________________________________________________________
623 AliTRDdataArray *AliTRDsegmentArray::GetDataArray(Int_t det) const
626 // Returns the data array for a given detector
629 return ((AliTRDdataArray *) AliTRDsegmentArray::At(det));
633 //_____________________________________________________________________________
634 AliTRDdataArray *AliTRDsegmentArray::GetDataArray(Int_t pla
639 // Returns the data array for a given detector
642 Int_t det = AliTRDgeometry::GetDetector(pla,cha,sec);
643 return GetDataArray(det);