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