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