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