]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDetector.cxx
User stepping methods added (E. Futo)
[u/mrichter/AliRoot.git] / STEER / AliDetector.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.22  2002/10/23 07:43:00  alibrary
19 Introducing some effective C++ suggestions
20
21 Revision 1.21  2002/10/22 15:02:15  alibrary
22 Introducing Riostream.h
23
24 Revision 1.20  2002/10/14 14:57:32  hristov
25 Merging the VirtualMC branch to the main development branch (HEAD)
26
27 Revision 1.16.8.3  2002/10/14 09:45:57  hristov
28 Updating VirtualMC to v3-09-02
29
30 Revision 1.19  2002/09/23 09:19:54  hristov
31 FirsTrackReference updated (M.Ivanov)
32
33 Revision 1.18  2002/08/26 13:51:17  hristov
34 Remaping of track references at the end of each primary particle (M.Ivanov)
35
36 Revision 1.17  2002/05/24 13:29:58  hristov
37 AliTrackReference added, AliDisplay modified
38
39 Revision 1.16  2001/10/04 15:30:56  hristov
40 Changes to accommodate the set of PHOS folders and tasks (Y.Schutz)
41
42 Revision 1.15  2001/07/27 13:03:13  hristov
43 Default Branch split level set to 99
44
45 Revision 1.14  2001/05/21 17:22:51  buncic
46 Fixed problem with missing AliConfig while reading galice.root
47
48 Revision 1.13  2001/05/16 14:57:22  alibrary
49 New files for folders and Stack
50
51 Revision 1.12  2001/03/12 17:47:03  hristov
52 Changes needed on Sun with CC 5.0
53
54 Revision 1.11  2001/01/26 19:58:46  hristov
55 Major upgrade of AliRoot code
56
57 Revision 1.10  2001/01/17 10:50:50  hristov
58 Corrections to destructors
59
60 Revision 1.9  2000/12/12 18:19:06  alibrary
61 Introduce consistency check when loading points
62
63 Revision 1.8  2000/11/30 07:12:48  alibrary
64 Introducing new Rndm and QA classes
65
66 Revision 1.7  2000/10/02 21:28:14  fca
67 Removal of useless dependecies via forward declarations
68
69 Revision 1.6  2000/07/12 08:56:25  fca
70 Coding convention correction and warning removal
71
72 Revision 1.5  1999/09/29 09:24:29  fca
73 Introduction of the Copyright and cvs Log
74
75 */
76
77 ///////////////////////////////////////////////////////////////////////////////
78 //                                                                           //
79 // Base class for ALICE modules. Both sensitive modules (detectors) and      //
80 // non-sensitive ones are described by this base class. This class           //
81 // supports the hit and digit trees produced by the simulation and also      //
82 // the objects produced by the reconstruction.                               //
83 //                                                                           //
84 // This class is also responsible for building the geometry of the           //
85 // detectors.                                                                //
86 //                                                                           //
87 //Begin_Html
88 /*
89 <img src="picts/AliDetectorClass.gif">
90 */
91 //End_Html
92 //                                                                           //
93 ///////////////////////////////////////////////////////////////////////////////
94
95 #include <assert.h>
96
97 #include <TBrowser.h>
98 #include <TFile.h>
99 #include <TFolder.h>
100 #include <TROOT.h>
101 #include <TTree.h>
102 #include <Riostream.h>
103
104 #include "AliConfig.h"
105 #include "AliDetector.h"
106 #include "AliHit.h"
107 #include "AliPoints.h"
108 #include "AliRun.h"
109 #include "AliTrackReference.h"
110
111
112 // Static variables for the hit iterator routines
113 static Int_t sMaxIterHit=0;
114 static Int_t sCurIterHit=0;
115
116
117 ClassImp(AliDetector)
118  
119 //_______________________________________________________________________
120 AliDetector::AliDetector():
121   fTimeGate(200.e-9),
122   fIshunt(0),
123   fNhits(0),
124   fNdigits(0),
125   fBufferSize(1600),
126   fHits(0),
127   fDigits(0),
128   fDigitsFile(0),
129   fPoints(0),
130   fTrackReferences(0),
131   fMaxIterTrackRef(0),
132   fCurrentIterTrackRef(0)
133 {
134   //
135   // Default constructor for the AliDetector class
136   //
137 }
138  
139 //_______________________________________________________________________
140 AliDetector::AliDetector(const AliDetector &det):
141   AliModule(det),
142   fTimeGate(200.e-9),
143   fIshunt(0),
144   fNhits(0),
145   fNdigits(0),
146   fBufferSize(1600),
147   fHits(0),
148   fDigits(0),
149   fDigitsFile(0),
150   fPoints(0),
151   fTrackReferences(0),
152   fMaxIterTrackRef(0),
153   fCurrentIterTrackRef(0)
154 {
155   det.Copy(*this);
156 }
157
158 //_____________________________________________________________________________
159 AliDetector::AliDetector(const char* name,const char *title):
160   AliModule(name,title),
161   fTimeGate(200.e-9),
162   fIshunt(0),
163   fNhits(0),
164   fNdigits(0),
165   fBufferSize(1600),
166   fHits(0),
167   fDigits(0),
168   fDigitsFile(0),
169   fPoints(0),
170   fTrackReferences(new TClonesArray("AliTrackReference", 100)),
171   fMaxIterTrackRef(0),
172   fCurrentIterTrackRef(0)
173 {
174   //
175   // Normal constructor invoked by all Detectors.
176   // Create the list for detector specific histograms
177   // Add this Detector to the global list of Detectors in Run.
178   //
179
180   fActive     = kTRUE;
181   AliConfig::Instance()->Add(this);
182
183 }
184  
185 //_______________________________________________________________________
186 AliDetector::~AliDetector()
187 {
188   //
189   // Destructor
190   //
191
192   // Delete space point structure
193   if (fPoints) {
194     fPoints->Delete();
195     delete fPoints;
196     fPoints     = 0;
197   }
198   // Delete digits structure
199   if (fDigits) {
200     fDigits->Delete();
201     delete fDigits;
202     fDigits     = 0;
203   }
204   if (fDigitsFile) delete [] fDigitsFile;
205 }
206
207 //_______________________________________________________________________
208 void AliDetector::Publish(const char *dir, void *address, const char *name)
209 {
210   //
211   // Register pointer to detector objects. 
212   // 
213   TFolder *topFolder = dynamic_cast<TFolder *>(gROOT->FindObjectAny("/Folders"));
214   if  (topFolder) { 
215     TFolder *folder = dynamic_cast<TFolder *>(topFolder->FindObjectAny(dir));
216     // TFolder *folder = dynamic_cast<TFolder *>(gROOT->FindObjectAny(dir));
217     if (!folder)  {
218       cerr << "Cannot register: Missing folder: " << dir << endl;
219     } else {
220       TFolder *subfolder = dynamic_cast<TFolder *>(folder->FindObjectAny(this->GetName())); 
221
222       if(!subfolder)
223          subfolder = folder->AddFolder(this->GetName(),this->GetTitle());
224       if (address) {
225         TObject **obj = static_cast<TObject **>(address);
226         if ((*obj)->InheritsFrom(TCollection::Class())) {
227            TCollection *collection = dynamic_cast<TCollection *>(*obj); 
228            if (name)
229              collection->SetName(name);
230         } 
231         subfolder->Add(*obj);
232       } 
233     }  
234   }
235 }
236
237 //_______________________________________________________________________
238 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, 
239                                        void* address, Int_t size,
240                                        const char *file)
241
242     return(MakeBranchInTree(tree,name,0,address,size,99,file));
243 }
244
245 //_______________________________________________________________________
246 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, 
247                                        const char *classname, 
248                                        void* address,Int_t size, 
249                                        Int_t splitlevel, const char *file)
250
251     //
252     // Makes branch in given tree and diverts them to a separate file
253     //  
254     if (GetDebug()>1)
255       printf("* MakeBranch * Making Branch %s \n",name);
256       
257     TDirectory *cwd = gDirectory;
258     TBranch *branch = 0;
259     
260     if (classname) {
261       branch = tree->Branch(name,classname,address,size,splitlevel);
262     } else {
263       branch = tree->Branch(name,address,size);
264     }
265        
266     if (file) {
267         char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
268         sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
269         branch->SetFile(outFile);
270         TIter next( branch->GetListOfBranches());
271         while ((branch=dynamic_cast<TBranch*>(next()))) {
272            branch->SetFile(outFile);
273         } 
274        delete outFile;
275         
276        cwd->cd();
277         
278        if (GetDebug()>1)
279            printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
280     }
281     const char *folder = 0;
282     TString folderName(name);  
283     
284     if (!strncmp(tree->GetName(),"TreeE",5)) folder = "RunMC/Event/Data";
285     if (!strncmp(tree->GetName(),"TreeK",5)) folder = "RunMC/Event/Data";
286     if (!strncmp(tree->GetName(),"TreeH",5)) {
287       folder     = "RunMC/Event/Data/Hits";
288       folderName = "Hits" ; 
289     }
290     if (!strncmp(tree->GetName(),"TreeTrackReferences",5)) {
291       folder     = "RunMC/Event/Data/TrackReferences";
292       folderName = "TrackReferences" ; 
293     }
294
295     if (!strncmp(tree->GetName(),"TreeD",5)) {
296       folder     = "Run/Event/Data";
297       folderName = "Digits" ; 
298     }
299     if (!strncmp(tree->GetName(),"TreeS",5)) {
300       folder     = "RunMC/Event/Data/SDigits";
301       folderName = "SDigits" ; 
302     }
303     if (!strncmp(tree->GetName(),"TreeR",5)) folder = "Run/Event/RecData";
304
305     if (folder) {
306       if (GetDebug())
307           printf("%15s: Publishing %s to %s\n",ClassName(),name,folder);
308       Publish(folder,address, folderName.Data());
309     }  
310     return branch;
311 }
312
313 //_______________________________________________________________________
314 void AliDetector::Browse(TBrowser *b)
315 {
316   //
317   // Insert Detector objects in the list of objects to be browsed
318   //
319   char name[64];
320   if( fHits == 0) return;
321   TObject *obj;
322   Int_t i, nobjects;
323   //
324   nobjects = fHits->GetEntries();
325   for (i=0;i<nobjects;i++) {
326     obj = fHits->At(i);
327     sprintf(name,"%s_%d",obj->GetName(),i);
328     b->Add(obj, &name[0]);
329   }
330 }
331
332 //_______________________________________________________________________
333 void AliDetector::Copy(AliDetector &) const
334 {
335   //
336   // Copy *this onto det -- not implemented
337   //
338   Fatal("Copy","Not implemented\n");
339 }
340
341 //_______________________________________________________________________
342 void AliDetector::FinishRun()
343 {
344   //
345   // Procedure called at the end of a run.
346   //
347 }
348
349 //_______________________________________________________________________
350 void AliDetector::RemapTrackReferencesIDs(Int_t *map)
351 {
352   // 
353   // Remapping track reference
354   // Called at finish primary
355   //
356   if (!fTrackReferences) return;
357   for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
358     AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
359     if (ref) {
360       Int_t newID = map[ref->GetTrack()];
361       if (newID>=0) ref->SetTrack(newID);
362       else ref->SetTrack(-1);
363       
364     }
365   }
366 }
367
368 //_______________________________________________________________________
369 AliHit* AliDetector::FirstHit(Int_t track)
370 {
371   //
372   // Initialise the hit iterator
373   // Return the address of the first hit for track
374   // If track>=0 the track is read from disk
375   // while if track<0 the first hit of the current
376   // track is returned
377   // 
378   if(track>=0) {
379     gAlice->ResetHits();
380     gAlice->TreeH()->GetEvent(track);
381   }
382   //
383   sMaxIterHit=fHits->GetEntriesFast();
384   sCurIterHit=0;
385   if(sMaxIterHit) return dynamic_cast<AliHit*>(fHits->UncheckedAt(0));
386   else            return 0;
387 }
388
389
390 //_______________________________________________________________________
391 AliTrackReference* AliDetector::FirstTrackReference(Int_t track)
392 {
393   //
394   // Initialise the hit iterator
395   // Return the address of the first hit for track
396   // If track>=0 the track is read from disk
397   // while if track<0 the first hit of the current
398   // track is returned
399   // 
400   if(track>=0) {
401     gAlice->ResetTrackReferences();
402     gAlice->TreeTR()->GetEvent(track);
403   }
404   //
405   fMaxIterTrackRef     = fTrackReferences->GetEntriesFast();
406   fCurrentIterTrackRef = 0;
407   if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
408   else            return 0;
409 }
410
411 //_______________________________________________________________________
412 AliHit* AliDetector::NextHit()
413 {
414   //
415   // Return the next hit for the current track
416   //
417   if(sMaxIterHit) {
418     if(++sCurIterHit<sMaxIterHit) 
419       return dynamic_cast<AliHit*>(fHits->UncheckedAt(sCurIterHit));
420     else        
421       return 0;
422   } else {
423     printf("* AliDetector::NextHit * Hit Iterator called without calling FistHit before\n");
424     return 0;
425   }
426 }
427
428 //_______________________________________________________________________
429 AliTrackReference* AliDetector::NextTrackReference()
430 {
431   //
432   // Return the next hit for the current track
433   //
434   if(fMaxIterTrackRef) {
435     if(++fCurrentIterTrackRef<fMaxIterTrackRef) 
436       return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
437     else        
438       return 0;
439   } else {
440     printf("* AliDetector::NextTrackReference * TrackReference  Iterator called without calling FistTrackReference before\n");
441     return 0;
442   }
443 }
444
445 //_______________________________________________________________________
446 void AliDetector::LoadPoints(Int_t)
447 {
448   //
449   // Store x, y, z of all hits in memory
450   //
451   if (fHits == 0) return;
452   //
453   Int_t nhits = fHits->GetEntriesFast();
454   if (nhits == 0) return;
455   Int_t tracks = gAlice->GetNtrack();
456   if (fPoints == 0) fPoints = new TObjArray(tracks);
457   AliHit *ahit;
458   //
459   Int_t *ntrk=new Int_t[tracks];
460   Int_t *limi=new Int_t[tracks];
461   Float_t **coor=new Float_t*[tracks];
462   for(Int_t i=0;i<tracks;i++) {
463     ntrk[i]=0;
464     coor[i]=0;
465     limi[i]=0;
466   }
467   //
468   AliPoints *points = 0;
469   Float_t *fp=0;
470   Int_t trk;
471   Int_t chunk=nhits/4+1;
472   //
473   // Loop over all the hits and store their position
474   for (Int_t hit=0;hit<nhits;hit++) {
475     ahit = dynamic_cast<AliHit*>(fHits->UncheckedAt(hit));
476     trk=ahit->GetTrack();
477     assert(trk<=tracks);
478     if(ntrk[trk]==limi[trk]) {
479       //
480       // Initialise a new track
481       fp=new Float_t[3*(limi[trk]+chunk)];
482       if(coor[trk]) {
483         memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
484         delete [] coor[trk];
485       }
486       limi[trk]+=chunk;
487       coor[trk] = fp;
488     } else {
489       fp = coor[trk];
490     }
491     fp[3*ntrk[trk]  ] = ahit->X();
492     fp[3*ntrk[trk]+1] = ahit->Y();
493     fp[3*ntrk[trk]+2] = ahit->Z();
494     ntrk[trk]++;
495   }
496   //
497   for(trk=0; trk<tracks; ++trk) {
498     if(ntrk[trk]) {
499       points = new AliPoints();
500       points->SetMarkerColor(GetMarkerColor());
501       points->SetMarkerSize(GetMarkerSize());
502       points->SetDetector(this);
503       points->SetParticle(trk);
504       points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
505       fPoints->AddAt(points,trk);
506       delete [] coor[trk];
507       coor[trk]=0;
508     }
509   }
510   delete [] coor;
511   delete [] ntrk;
512   delete [] limi;
513 }
514
515 //_______________________________________________________________________
516 void AliDetector::MakeBranch(Option_t *option, const char *file)
517 {
518   //
519   // Create a new branch in the current Root Tree
520   // The branch of fHits is automatically split
521   //
522  
523   char branchname[10];
524   sprintf(branchname,"%s",GetName());
525   //
526   // Get the pointer to the header
527   const char *cH = strstr(option,"H");
528   //
529   if (fHits && gAlice->TreeH() && cH) {
530     MakeBranchInTree(gAlice->TreeH(), 
531                      branchname, &fHits, fBufferSize, file) ;              
532   }     
533   
534   const char *cD = strstr(option,"D");
535
536   if (cD) {
537     if (file) {
538        fDigitsFile = new char[strlen (file)];
539        strcpy(fDigitsFile,file);
540     }
541   }
542 }
543 //_______________________________________________________________________
544 void AliDetector::MakeBranchTR(Option_t *option, const char *file)
545 {
546   //
547   // Create a new branch in the current Root Tree
548   // The branch of fHits is automatically split
549   //
550  
551   char branchname[10];
552   sprintf(branchname,"%s",GetName());
553   //
554   // Get the pointer to the header
555   const char *cTR = strstr(option,"T");
556   //
557   if (fTrackReferences && gAlice->TreeTR() && cTR) {
558     MakeBranchInTree(gAlice->TreeTR(), 
559                      branchname, &fTrackReferences, fBufferSize, file) ;              
560   }       
561 }
562
563 //_______________________________________________________________________
564 void AliDetector::ResetDigits()
565 {
566   //
567   // Reset number of digits and the digits array
568   //
569   fNdigits   = 0;
570   if (fDigits)   fDigits->Clear();
571 }
572
573 //_______________________________________________________________________
574 void AliDetector::ResetHits()
575 {
576   //
577   // Reset number of hits and the hits array
578   //
579   fNhits   = 0;
580   if (fHits)   fHits->Clear();
581 }
582
583 //_______________________________________________________________________
584 void AliDetector::ResetTrackReferences()
585 {
586   //
587   // Reset number of hits and the hits array
588   //
589   fMaxIterTrackRef   = 0;
590   if (fTrackReferences)   fTrackReferences->Clear();
591 }
592
593 //_______________________________________________________________________
594 void AliDetector::ResetPoints()
595 {
596   //
597   // Reset array of points
598   //
599   if (fPoints) {
600     fPoints->Delete();
601     delete fPoints;
602     fPoints = 0;
603   }
604 }
605
606 //_______________________________________________________________________
607 void AliDetector::SetTreeAddress()
608 {
609   //
610   // Set branch address for the Hits and Digits Trees
611   //
612   TBranch *branch;
613   char branchname[20];
614   sprintf(branchname,"%s",GetName());
615   //
616   // Branch address for hit tree
617   TTree *treeH = gAlice->TreeH();
618   if (treeH && fHits) {
619     branch = treeH->GetBranch(branchname);
620     if (branch) branch->SetAddress(&fHits);
621   }
622   //
623   // Branch address for digit tree
624   TTree *treeD = gAlice->TreeD();
625   if (treeD && fDigits) {
626     branch = treeD->GetBranch(branchname);
627     if (branch) branch->SetAddress(&fDigits);
628   }
629
630   // Branch address for tr  tree
631   TTree *treeTR = gAlice->TreeTR();
632   if (treeTR && fTrackReferences) {
633     branch = treeTR->GetBranch(branchname);
634     if (branch) branch->SetAddress(&fTrackReferences);
635   }
636 }
637
638