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