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