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