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