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