]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDetector.cxx
Code for simulation, sdigitization and digitization moved from macros to compiled...
[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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // Base class for ALICE modules. Both sensitive modules (detectors) and      //
21 // non-sensitive ones are described by this base class. This class           //
22 // supports the hit and digit trees produced by the simulation and also      //
23 // the objects produced by the reconstruction.                               //
24 //                                                                           //
25 // This class is also responsible for building the geometry of the           //
26 // detectors.                                                                //
27 //                                                                           //
28 //Begin_Html
29 /*
30 <img src="picts/AliDetectorClass.gif">
31 */
32 //End_Html
33 //                                                                           //
34 ///////////////////////////////////////////////////////////////////////////////
35
36 #include <assert.h>
37
38 #include <TBrowser.h>
39 #include <TFile.h>
40 #include <TFolder.h>
41 #include <TROOT.h>
42 #include <TTree.h>
43 #include <Riostream.h>
44
45 #include "AliConfig.h"
46 #include "AliDetector.h"
47 #include "AliHit.h"
48 #include "AliPoints.h"
49 #include "AliLoader.h"
50 #include "AliRun.h"
51 #include "AliMC.h"
52
53
54 // Static variables for the hit iterator routines
55 static Int_t sMaxIterHit=0;
56 static Int_t sCurIterHit=0;
57
58
59 ClassImp(AliDetector)
60  
61 //_______________________________________________________________________
62 AliDetector::AliDetector():
63   fTimeGate(200.e-9),
64   fIshunt(0),
65   fNhits(0),
66   fNdigits(0),
67   fBufferSize(1600),
68   fHits(0),
69   fDigits(0),
70   fPoints(0),
71   fLoader(0x0)
72 {
73   //
74   // Default constructor for the AliDetector class
75   //
76 }
77  
78 //_______________________________________________________________________
79 AliDetector::AliDetector(const AliDetector &det):
80   AliModule(det),
81   fTimeGate(200.e-9),
82   fIshunt(0),
83   fNhits(0),
84   fNdigits(0),
85   fBufferSize(1600),
86   fHits(0),
87   fDigits(0),
88   fPoints(0),
89   fLoader(0x0)
90 {
91   det.Copy(*this);
92 }
93
94 //_____________________________________________________________________________
95 AliDetector::AliDetector(const char* name,const char *title):
96   AliModule(name,title),
97   fTimeGate(200.e-9),
98   fIshunt(0),
99   fNhits(0),
100   fNdigits(0),
101   fBufferSize(1600),
102   fHits(0),
103   fDigits(0),
104   fPoints(0),
105   fLoader(0x0)
106 {
107   //
108   // Normal constructor invoked by all Detectors.
109   // Create the list for detector specific histograms
110   // Add this Detector to the global list of Detectors in Run.
111   //
112
113   fActive     = kTRUE;
114   AliConfig::Instance()->Add(this);
115
116 }
117  
118 //_______________________________________________________________________
119 AliDetector::~AliDetector()
120 {
121   //
122   // Destructor
123   //
124
125   // Delete space point structure
126   if (fPoints) {
127     fPoints->Delete();
128     delete fPoints;
129     fPoints     = 0;
130   }
131   // Delete digits structure
132   if (fDigits) {
133     fDigits->Delete();
134     delete fDigits;
135     fDigits     = 0;
136   }
137
138   if (fLoader)
139    {
140     fLoader->GetModulesFolder()->Remove(this);
141    }
142
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   MayNotUse("Publish");
153 }
154
155 //_______________________________________________________________________
156 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, 
157                                        void* address, Int_t size,
158                                        const char *file)
159
160     return(MakeBranchInTree(tree,name,0,address,size,99,file));
161 }
162
163 //_______________________________________________________________________
164 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, 
165                                        const char *classname, 
166                                        void* address,Int_t size, 
167                                        Int_t splitlevel, const char */*file*/)
168
169 //
170 // Makes branch in given tree and diverts them to a separate file
171 // 
172 //
173 //
174 // if (GetDebug()>1)
175  if(GetDebug()) Info("MakeBranch","Making Branch %s",name);
176  if (tree == 0x0) 
177   {
178    Error("MakeBranch","Making Branch %s Tree is NULL",name);
179    return 0x0;
180   }
181  TBranch *branch = tree->GetBranch(name);
182  if (branch) 
183   {  
184     if(GetDebug()) Info("MakeBranch","Branch %s is already in tree.",name);
185     return branch;
186   }
187     
188  if (classname) 
189   {
190     branch = tree->Branch(name,classname,address,size,splitlevel);
191   } 
192  else 
193   {
194     branch = tree->Branch(name,address,size);
195   }
196  if(GetDebug()) Info("MakeBranch","Branch %s returning branch %#x",name,branch);
197  return branch;
198 }
199
200 //_______________________________________________________________________
201 void AliDetector::Browse(TBrowser *b)
202 {
203   //
204   // Insert Detector objects in the list of objects to be browsed
205   //
206   char name[64];
207   if( fHits == 0) return;
208   TObject *obj;
209   Int_t i, nobjects;
210   //
211   nobjects = fHits->GetEntries();
212   for (i=0;i<nobjects;i++) {
213     obj = fHits->At(i);
214     sprintf(name,"%s_%d",obj->GetName(),i);
215     b->Add(obj, &name[0]);
216   }
217 }
218
219 //_______________________________________________________________________
220 void AliDetector::Copy(AliDetector &) const
221 {
222   //
223   // Copy *this onto det -- not implemented
224   //
225   Fatal("Copy","Not implemented\n");
226 }
227
228 //_______________________________________________________________________
229 void AliDetector::FinishRun()
230 {
231   //
232   // Procedure called at the end of a run.
233   //
234 }
235
236 //_______________________________________________________________________
237 AliHit* AliDetector::FirstHit(Int_t track)
238 {
239   //
240   // Initialise the hit iterator
241   // Return the address of the first hit for track
242   // If track>=0 the track is read from disk
243   // while if track<0 the first hit of the current
244   // track is returned
245   // 
246   if(track>=0) {
247     gAlice->ResetHits(); //stupid = if N detector this method is called N times
248     TreeH()->GetEvent(track); //skowron
249   }
250   //
251   sMaxIterHit=fHits->GetEntriesFast();
252   sCurIterHit=0;
253   if(sMaxIterHit) return dynamic_cast<AliHit*>(fHits->UncheckedAt(0));
254   else            return 0;
255 }
256
257 //_______________________________________________________________________
258 AliHit* AliDetector::NextHit()
259 {
260   //
261   // Return the next hit for the current track
262   //
263   if(sMaxIterHit) {
264     if(++sCurIterHit<sMaxIterHit) 
265       return dynamic_cast<AliHit*>(fHits->UncheckedAt(sCurIterHit));
266     else        
267       return 0;
268   } else {
269     printf("* AliDetector::NextHit * Hit Iterator called without calling FistHit before\n");
270     return 0;
271   }
272 }
273
274 //_______________________________________________________________________
275 void AliDetector::LoadPoints(Int_t)
276 {
277   //
278   // Store x, y, z of all hits in memory
279   //
280   if (fHits == 0) 
281    {
282     Error("LoadPoints","fHits == 0. Name is %s",GetName());
283     return;
284    }
285   //
286   Int_t nhits = fHits->GetEntriesFast();
287   if (nhits == 0) 
288    {
289 //    Error("LoadPoints","nhits == 0. Name is %s",GetName());
290     return;
291    }
292   Int_t tracks = gAlice->GetMCApp()->GetNtrack();
293   if (fPoints == 0) fPoints = new TObjArray(tracks);
294   AliHit *ahit;
295   //
296   Int_t *ntrk=new Int_t[tracks];
297   Int_t *limi=new Int_t[tracks];
298   Float_t **coor=new Float_t*[tracks];
299   for(Int_t i=0;i<tracks;i++) {
300     ntrk[i]=0;
301     coor[i]=0;
302     limi[i]=0;
303   }
304   //
305   AliPoints *points = 0;
306   Float_t *fp=0;
307   Int_t trk;
308   Int_t chunk=nhits/4+1;
309   //
310   // Loop over all the hits and store their position
311   for (Int_t hit=0;hit<nhits;hit++) {
312     ahit = dynamic_cast<AliHit*>(fHits->UncheckedAt(hit));
313     trk=ahit->GetTrack();
314     assert(trk<=tracks);
315     if(ntrk[trk]==limi[trk])
316      {
317       //
318       // Initialise a new track
319       fp=new Float_t[3*(limi[trk]+chunk)];
320       if(coor[trk]) 
321        {
322           memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
323           delete [] coor[trk];
324        }
325       limi[trk]+=chunk;
326       coor[trk] = fp;
327      } 
328     else 
329      {
330       fp = coor[trk];
331      }
332     fp[3*ntrk[trk]  ] = ahit->X();
333     fp[3*ntrk[trk]+1] = ahit->Y();
334     fp[3*ntrk[trk]+2] = ahit->Z();
335     ntrk[trk]++;
336   }
337   //
338   for(trk=0; trk<tracks; ++trk) {
339     if(ntrk[trk]) {
340       points = new AliPoints();
341       points->SetMarkerColor(GetMarkerColor());
342       points->SetMarkerSize(GetMarkerSize());
343       points->SetDetector(this);
344       points->SetParticle(trk);
345       points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
346       fPoints->AddAt(points,trk);
347       delete [] coor[trk];
348       coor[trk]=0;
349     }
350   }
351   delete [] coor;
352   delete [] ntrk;
353   delete [] limi;
354 }
355
356 //_______________________________________________________________________
357 void AliDetector::MakeBranch(Option_t *option)
358 {
359   //
360   // Create a new branch for this detector in its treeH
361   //
362
363   if(GetDebug()) Info("MakeBranch"," for %s",GetName());
364   const char *cH = strstr(option,"H");
365
366   if (fHits && TreeH() && cH) 
367    {
368      MakeBranchInTree(TreeH(), GetName(), &fHits, fBufferSize, 0);
369    }    
370 }
371
372 //_______________________________________________________________________
373 void AliDetector::ResetDigits()
374 {
375   //
376   // Reset number of digits and the digits array
377   //
378   fNdigits   = 0;
379   if (fDigits) fDigits->Clear();
380 }
381
382 //_______________________________________________________________________
383 void AliDetector::ResetHits()
384 {
385   //
386   // Reset number of hits and the hits array
387   //
388   fNhits   = 0;
389   if (fHits) fHits->Clear();
390 }
391
392 //_______________________________________________________________________
393 void AliDetector::ResetPoints()
394 {
395   //
396   // Reset array of points
397   //
398   if (fPoints) {
399     fPoints->Delete();
400     delete fPoints;
401     fPoints = 0;
402   }
403 }
404
405 //_______________________________________________________________________
406 void AliDetector::SetTreeAddress()
407 {
408   //
409   // Set branch address for the Hits and Digits Trees
410   //
411   TBranch *branch;
412   //
413   // Branch address for hit tree
414   
415   TTree *tree = TreeH();
416   if (tree && fHits) {
417     branch = tree->GetBranch(GetName());
418     if (branch) 
419      {
420        if(GetDebug()) Info("SetTreeAddress","(%s) Setting for Hits",GetName());
421        branch->SetAddress(&fHits);
422      }
423     else
424      { //can be invoked before branch creation
425        if(GetDebug()) Warning("SetTreeAddress","(%s) Failed for Hits. Can not find branch in tree.",GetName());
426      }
427   }
428   
429   //
430   // Branch address for digit tree
431   TTree *treeD = fLoader->TreeD();
432   if (treeD && fDigits) {
433     branch = treeD->GetBranch(GetName());
434     if (branch) branch->SetAddress(&fDigits);
435   }
436   
437   AliModule::SetTreeAddress();
438 }
439
440 //_______________________________________________________________________
441 void AliDetector::MakeTree(Option_t *option)
442  {
443  //makes a tree (container) for the data defined in option
444  //"H" - hits
445  //"D" - digits
446  //"S" - summable digits
447  //"R" - recontructed points and tracks
448  
449     AliLoader* loader = GetLoader();
450     if (loader == 0x0)
451      {
452        Error("MakeTree","Can not get loader for %s",GetName());
453        return;
454      }
455     loader->MakeTree(option); //delegate this job to getter
456  }
457
458 //_______________________________________________________________________
459 AliLoader* AliDetector::MakeLoader(const char* topfoldername)
460
461 //builds standard getter (AliLoader type)
462 //if detector wants to use castomized getter, it must overload this method
463
464  if (GetDebug())
465    Info("MakeLoader",
466         "Creating standard getter for detector %s. Top folder is %s.",
467          GetName(),topfoldername);
468      
469  fLoader = new AliLoader(GetName(),topfoldername);
470  return fLoader;
471 }
472
473 //_______________________________________________________________________
474 TTree* AliDetector::TreeH()
475 {
476 //Get the hits container from the folder
477   if (GetLoader() == 0x0) 
478     {
479     //sunstitude this with make getter when we can obtain the event folder name 
480      Error("TreeH","Can not get the getter");
481      return 0x0;
482     }
483  
484   TTree* tree = (TTree*)GetLoader()->TreeH();
485   return tree;
486 }
487