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