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