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