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 Revision 1.18 2002/08/26 13:51:17 hristov
19 Remaping of track references at the end of each primary particle (M.Ivanov)
21 Revision 1.17 2002/05/24 13:29:58 hristov
22 AliTrackReference added, AliDisplay modified
24 Revision 1.16 2001/10/04 15:30:56 hristov
25 Changes to accommodate the set of PHOS folders and tasks (Y.Schutz)
27 Revision 1.15 2001/07/27 13:03:13 hristov
28 Default Branch split level set to 99
30 Revision 1.14 2001/05/21 17:22:51 buncic
31 Fixed problem with missing AliConfig while reading galice.root
33 Revision 1.13 2001/05/16 14:57:22 alibrary
34 New files for folders and Stack
36 Revision 1.12 2001/03/12 17:47:03 hristov
37 Changes needed on Sun with CC 5.0
39 Revision 1.11 2001/01/26 19:58:46 hristov
40 Major upgrade of AliRoot code
42 Revision 1.10 2001/01/17 10:50:50 hristov
43 Corrections to destructors
45 Revision 1.9 2000/12/12 18:19:06 alibrary
46 Introduce consistency check when loading points
48 Revision 1.8 2000/11/30 07:12:48 alibrary
49 Introducing new Rndm and QA classes
51 Revision 1.7 2000/10/02 21:28:14 fca
52 Removal of useless dependecies via forward declarations
54 Revision 1.6 2000/07/12 08:56:25 fca
55 Coding convention correction and warning removal
57 Revision 1.5 1999/09/29 09:24:29 fca
58 Introduction of the Copyright and cvs Log
62 ///////////////////////////////////////////////////////////////////////////////
64 // Base class for ALICE modules. Both sensitive modules (detectors) and //
65 // non-sensitive ones are described by this base class. This class //
66 // supports the hit and digit trees produced by the simulation and also //
67 // the objects produced by the reconstruction. //
69 // This class is also responsible for building the geometry of the //
74 <img src="picts/AliDetectorClass.gif">
78 ///////////////////////////////////////////////////////////////////////////////
89 #include "AliConfig.h"
90 #include "AliDetector.h"
93 #include "AliPoints.h"
94 #include "AliTrackReference.h"
97 // Static variables for the hit iterator routines
98 static Int_t sMaxIterHit=0;
99 static Int_t sCurIterHit=0;
102 ClassImp(AliDetector)
104 //_____________________________________________________________________________
105 AliDetector::AliDetector()
108 // Default constructor for the AliDetector class
121 //_____________________________________________________________________________
122 AliDetector::AliDetector(const char* name,const char *title):AliModule(name,title)
125 // Normal constructor invoked by all Detectors.
126 // Create the list for detector specific histograms
127 // Add this Detector to the global list of Detectors in Run.
141 AliConfig::Instance()->Add(this);
143 fTrackReferences = new TClonesArray("AliTrackReference", 100);
144 //if detector to want to create another track reference - than let's be free
148 //_____________________________________________________________________________
149 AliDetector::~AliDetector()
157 // Delete space point structure
163 // Delete digits structure
169 if (fDigitsFile) delete [] fDigitsFile;
172 //_____________________________________________________________________________
173 void AliDetector::Publish(const char *dir, void *address, const char *name)
176 // Register pointer to detector objects.
178 TFolder *topFolder = (TFolder *)gROOT->FindObjectAny("/Folders");
180 TFolder *folder = (TFolder *)topFolder->FindObjectAny(dir);
181 // TFolder *folder = (TFolder *)gROOT->FindObjectAny(dir);
183 cerr << "Cannot register: Missing folder: " << dir << endl;
185 TFolder *subfolder = (TFolder *) folder->FindObjectAny(this->GetName());
188 subfolder = folder->AddFolder(this->GetName(),this->GetTitle());
190 TObject **obj = (TObject **) address;
191 if ((*obj)->InheritsFrom(TCollection::Class())) {
192 TCollection *collection = (TCollection *) (*obj);
194 collection->SetName(name);
196 subfolder->Add(*obj);
202 //_____________________________________________________________________________
203 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, void* address, Int_t size,const char *file)
205 return(MakeBranchInTree(tree,name,0,address,size,99,file));
208 //_____________________________________________________________________________
209 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, const char *classname, void* address,Int_t size, Int_t splitlevel, const char *file)
212 // Makes branch in given tree and diverts them to a separate file
215 printf("* MakeBranch * Making Branch %s \n",name);
217 TDirectory *cwd = gDirectory;
221 branch = tree->Branch(name,classname,address,size,splitlevel);
223 branch = tree->Branch(name,address,size);
227 char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
228 sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
229 branch->SetFile(outFile);
230 TIter next( branch->GetListOfBranches());
231 while ((branch=(TBranch*)next())) {
232 branch->SetFile(outFile);
239 printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
242 TString folderName(name);
244 if (!strncmp(tree->GetName(),"TreeE",5)) folder = "RunMC/Event/Data";
245 if (!strncmp(tree->GetName(),"TreeK",5)) folder = "RunMC/Event/Data";
246 if (!strncmp(tree->GetName(),"TreeH",5)) {
247 folder = "RunMC/Event/Data/Hits";
248 folderName = "Hits" ;
250 if (!strncmp(tree->GetName(),"TreeTrackReferences",5)) {
251 folder = "RunMC/Event/Data/TrackReferences";
252 folderName = "TrackReferences" ;
255 if (!strncmp(tree->GetName(),"TreeD",5)) {
256 folder = "Run/Event/Data";
257 folderName = "Digits" ;
259 if (!strncmp(tree->GetName(),"TreeS",5)) {
260 folder = "RunMC/Event/Data/SDigits";
261 folderName = "SDigits" ;
263 if (!strncmp(tree->GetName(),"TreeR",5)) folder = "Run/Event/RecData";
267 printf("%15s: Publishing %s to %s\n",ClassName(),name,folder);
268 Publish(folder,address, folderName.Data());
273 //_____________________________________________________________________________
274 void AliDetector::Browse(TBrowser *b)
277 // Insert Detector objects in the list of objects to be browsed
280 if( fHits == 0) return;
284 nobjects = fHits->GetEntries();
285 for (i=0;i<nobjects;i++) {
287 sprintf(name,"%s_%d",obj->GetName(),i);
288 b->Add(obj, &name[0]);
292 //_____________________________________________________________________________
293 void AliDetector::Copy(AliDetector &det) const
296 // Copy *this onto det -- not implemented
298 Fatal("Copy","Not implemented~\n");
301 //_____________________________________________________________________________
302 void AliDetector::FinishRun()
305 // Procedure called at the end of a run.
311 void AliDetector::RemapTrackReferencesIDs(Int_t *map)
314 //remaping track reference
315 //called at finish primary
316 if (!fTrackReferences) return;
317 for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
318 AliTrackReference * ref = (AliTrackReference*) fTrackReferences->UncheckedAt(i);
320 Int_t newID = map[ref->GetTrack()];
321 if (newID>=0) ref->SetTrack(newID);
322 else ref->SetTrack(-1);
328 //_____________________________________________________________________________
329 AliHit* AliDetector::FirstHit(Int_t track)
332 // Initialise the hit iterator
333 // Return the address of the first hit for track
334 // If track>=0 the track is read from disk
335 // while if track<0 the first hit of the current
340 gAlice->TreeH()->GetEvent(track);
343 sMaxIterHit=fHits->GetEntriesFast();
345 if(sMaxIterHit) return (AliHit*) fHits->UncheckedAt(0);
350 //_____________________________________________________________________________
351 AliTrackReference* AliDetector::FirstTrackReference(Int_t track)
354 // Initialise the hit iterator
355 // Return the address of the first hit for track
356 // If track>=0 the track is read from disk
357 // while if track<0 the first hit of the current
361 gAlice->ResetTrackReferences();
362 gAlice->TreeTR()->GetEvent(track);
365 fMaxIterTrackRef = fTrackReferences->GetEntriesFast();
366 fCurrentIterTrackRef = 0;
367 if(fMaxIterTrackRef) return (AliTrackReference*) fTrackReferences->UncheckedAt(0);
373 //_____________________________________________________________________________
374 AliHit* AliDetector::NextHit()
377 // Return the next hit for the current track
380 if(++sCurIterHit<sMaxIterHit)
381 return (AliHit*) fHits->UncheckedAt(sCurIterHit);
385 printf("* AliDetector::NextHit * Hit Iterator called without calling FistHit before\n");
389 //_____________________________________________________________________________
390 AliTrackReference* AliDetector::NextTrackReference()
393 // Return the next hit for the current track
395 if(fMaxIterTrackRef) {
396 if(++fCurrentIterTrackRef<fMaxIterTrackRef)
397 return (AliTrackReference*) fTrackReferences->UncheckedAt(fCurrentIterTrackRef);
401 printf("* AliDetector::NextTrackReference * TrackReference Iterator called without calling FistTrackReference before\n");
406 //_____________________________________________________________________________
407 void AliDetector::LoadPoints(Int_t)
410 // Store x, y, z of all hits in memory
412 if (fHits == 0) return;
414 Int_t nhits = fHits->GetEntriesFast();
415 if (nhits == 0) return;
416 Int_t tracks = gAlice->GetNtrack();
417 if (fPoints == 0) fPoints = new TObjArray(tracks);
420 Int_t *ntrk=new Int_t[tracks];
421 Int_t *limi=new Int_t[tracks];
422 Float_t **coor=new Float_t*[tracks];
423 for(Int_t i=0;i<tracks;i++) {
429 AliPoints *points = 0;
432 Int_t chunk=nhits/4+1;
434 // Loop over all the hits and store their position
435 for (Int_t hit=0;hit<nhits;hit++) {
436 ahit = (AliHit*)fHits->UncheckedAt(hit);
437 trk=ahit->GetTrack();
439 if(ntrk[trk]==limi[trk]) {
441 // Initialise a new track
442 fp=new Float_t[3*(limi[trk]+chunk)];
444 memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
452 fp[3*ntrk[trk] ] = ahit->X();
453 fp[3*ntrk[trk]+1] = ahit->Y();
454 fp[3*ntrk[trk]+2] = ahit->Z();
458 for(trk=0; trk<tracks; ++trk) {
460 points = new AliPoints();
461 points->SetMarkerColor(GetMarkerColor());
462 points->SetMarkerSize(GetMarkerSize());
463 points->SetDetector(this);
464 points->SetParticle(trk);
465 points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
466 fPoints->AddAt(points,trk);
476 //_____________________________________________________________________________
477 void AliDetector::MakeBranch(Option_t *option, const char *file)
480 // Create a new branch in the current Root Tree
481 // The branch of fHits is automatically split
485 sprintf(branchname,"%s",GetName());
487 // Get the pointer to the header
488 const char *cH = strstr(option,"H");
490 if (fHits && gAlice->TreeH() && cH) {
491 MakeBranchInTree(gAlice->TreeH(),
492 branchname, &fHits, fBufferSize, file) ;
495 const char *cD = strstr(option,"D");
499 fDigitsFile = new char[strlen (file)];
500 strcpy(fDigitsFile,file);
504 //_____________________________________________________________________________
505 void AliDetector::MakeBranchTR(Option_t *option, const char *file)
508 // Create a new branch in the current Root Tree
509 // The branch of fHits is automatically split
513 sprintf(branchname,"%s",GetName());
515 // Get the pointer to the header
516 const char *cTR = strstr(option,"T");
518 if (fTrackReferences && gAlice->TreeTR() && cTR) {
519 MakeBranchInTree(gAlice->TreeTR(),
520 branchname, &fTrackReferences, fBufferSize, file) ;
524 //_____________________________________________________________________________
525 void AliDetector::ResetDigits()
528 // Reset number of digits and the digits array
531 if (fDigits) fDigits->Clear();
534 //_____________________________________________________________________________
535 void AliDetector::ResetHits()
538 // Reset number of hits and the hits array
541 if (fHits) fHits->Clear();
546 //_____________________________________________________________________________
547 void AliDetector::ResetTrackReferences()
550 // Reset number of hits and the hits array
552 fMaxIterTrackRef = 0;
553 if (fTrackReferences) fTrackReferences->Clear();
558 //_____________________________________________________________________________
559 void AliDetector::ResetPoints()
562 // Reset array of points
571 //_____________________________________________________________________________
572 void AliDetector::SetTreeAddress()
575 // Set branch address for the Hits and Digits Trees
579 sprintf(branchname,"%s",GetName());
581 // Branch address for hit tree
582 TTree *treeH = gAlice->TreeH();
583 if (treeH && fHits) {
584 branch = treeH->GetBranch(branchname);
585 if (branch) branch->SetAddress(&fHits);
588 // Branch address for digit tree
589 TTree *treeD = gAlice->TreeD();
590 if (treeD && fDigits) {
591 branch = treeD->GetBranch(branchname);
592 if (branch) branch->SetAddress(&fDigits);
595 // Branch address for tr tree
596 TTree *treeTR = gAlice->TreeTR();
597 if (treeTR && fTrackReferences) {
598 branch = treeTR->GetBranch(branchname);
599 if (branch) branch->SetAddress(&fTrackReferences);