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