]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDetector.cxx
Updated class version
[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 <assert.h>
37
38 #include <Riostream.h>
39 #include <TBrowser.h>
40 #include <TFile.h>
41 #include <TFolder.h>
42 #include <TROOT.h>
43 #include <TTree.h>
44
45 #include "AliConfig.h"
46 #include "AliDetector.h"
47 #include "AliHit.h"
48 #include "AliPoints.h"
49 #include "AliRun.h"
50 #include "AliTrackReference.h"
51
52
53 // Static variables for the hit iterator routines
54 static Int_t sMaxIterHit=0;
55 static Int_t sCurIterHit=0;
56
57
58 ClassImp(AliDetector)
59  
60 //_______________________________________________________________________
61 AliDetector::AliDetector():
62   fTimeGate(200.e-9),
63   fIshunt(0),
64   fNhits(0),
65   fNdigits(0),
66   fBufferSize(1600),
67   fHits(0),
68   fDigits(0),
69   fDigitsFile(0),
70   fPoints(0),
71   fTrackReferences(0),
72   fMaxIterTrackRef(0),
73   fCurrentIterTrackRef(0)
74 {
75   //
76   // Default constructor for the AliDetector class
77   //
78 }
79  
80 //_______________________________________________________________________
81 AliDetector::AliDetector(const AliDetector &det):
82   AliModule(det),
83   fTimeGate(200.e-9),
84   fIshunt(0),
85   fNhits(0),
86   fNdigits(0),
87   fBufferSize(1600),
88   fHits(0),
89   fDigits(0),
90   fDigitsFile(0),
91   fPoints(0),
92   fTrackReferences(0),
93   fMaxIterTrackRef(0),
94   fCurrentIterTrackRef(0)
95 {
96   det.Copy(*this);
97 }
98
99 //_____________________________________________________________________________
100 AliDetector::AliDetector(const char* name,const char *title):
101   AliModule(name,title),
102   fTimeGate(200.e-9),
103   fIshunt(0),
104   fNhits(0),
105   fNdigits(0),
106   fBufferSize(1600),
107   fHits(0),
108   fDigits(0),
109   fDigitsFile(0),
110   fPoints(0),
111   fTrackReferences(new TClonesArray("AliTrackReference", 100)),
112   fMaxIterTrackRef(0),
113   fCurrentIterTrackRef(0)
114 {
115   //
116   // Normal constructor invoked by all Detectors.
117   // Create the list for detector specific histograms
118   // Add this Detector to the global list of Detectors in Run.
119   //
120
121   fActive     = kTRUE;
122   AliConfig::Instance()->Add(this);
123
124 }
125  
126 //_______________________________________________________________________
127 AliDetector::~AliDetector()
128 {
129   //
130   // Destructor
131   //
132
133   // Delete space point structure
134   if (fPoints) {
135     fPoints->Delete();
136     delete fPoints;
137     fPoints     = 0;
138   }
139   // Delete digits structure
140   if (fDigits) {
141     fDigits->Delete();
142     delete fDigits;
143     fDigits     = 0;
144   }
145   if (fDigitsFile) delete [] fDigitsFile;
146 }
147
148 //_______________________________________________________________________
149 void AliDetector::Publish(const char *dir, void *address, const char *name)
150 {
151   //
152   // Register pointer to detector objects. 
153   // 
154   TFolder *topFolder = dynamic_cast<TFolder *>(gROOT->FindObjectAny("/Folders"));
155   if  (topFolder) { 
156     TFolder *folder = dynamic_cast<TFolder *>(topFolder->FindObjectAny(dir));
157     // TFolder *folder = dynamic_cast<TFolder *>(gROOT->FindObjectAny(dir));
158     if (!folder)  {
159       cerr << "Cannot register: Missing folder: " << dir << endl;
160     } else {
161       TFolder *subfolder = dynamic_cast<TFolder *>(folder->FindObjectAny(this->GetName())); 
162
163       if(!subfolder)
164          subfolder = folder->AddFolder(this->GetName(),this->GetTitle());
165       if (address) {
166         TObject **obj = static_cast<TObject **>(address);
167         if ((*obj)->InheritsFrom(TCollection::Class())) {
168            TCollection *collection = dynamic_cast<TCollection *>(*obj); 
169            if (name)
170              collection->SetName(name);
171         } 
172         subfolder->Add(*obj);
173       } 
174     }  
175   }
176 }
177
178 //_______________________________________________________________________
179 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, 
180                                        void* address, Int_t size,
181                                        const char *file)
182
183     return(MakeBranchInTree(tree,name,0,address,size,99,file));
184 }
185
186 //_______________________________________________________________________
187 TBranch* AliDetector::MakeBranchInTree(TTree *tree, const char* name, 
188                                        const char *classname, 
189                                        void* address,Int_t size, 
190                                        Int_t splitlevel, const char *file)
191
192     //
193     // Makes branch in given tree and diverts them to a separate file
194     //  
195     if (GetDebug()>1)
196       printf("* MakeBranch * Making Branch %s \n",name);
197       
198     TDirectory *cwd = gDirectory;
199     TBranch *branch = 0;
200     
201     if (classname) {
202       branch = tree->Branch(name,classname,address,size,splitlevel);
203     } else {
204       branch = tree->Branch(name,address,size);
205     }
206        
207     if (file) {
208         char * outFile = new char[strlen(gAlice->GetBaseFile())+strlen(file)+2];
209         sprintf(outFile,"%s/%s",gAlice->GetBaseFile(),file);
210         branch->SetFile(outFile);
211         TIter next( branch->GetListOfBranches());
212         while ((branch=dynamic_cast<TBranch*>(next()))) {
213            branch->SetFile(outFile);
214         } 
215        delete outFile;
216         
217        cwd->cd();
218         
219        if (GetDebug()>1)
220            printf("* MakeBranch * Diverting Branch %s to file %s\n",name,file);
221     }
222     const char *folder = 0;
223     TString folderName(name);  
224     
225     if (!strncmp(tree->GetName(),"TreeE",5)) folder = "RunMC/Event/Data";
226     if (!strncmp(tree->GetName(),"TreeK",5)) folder = "RunMC/Event/Data";
227     if (!strncmp(tree->GetName(),"TreeH",5)) {
228       folder     = "RunMC/Event/Data/Hits";
229       folderName = "Hits" ; 
230     }
231     if (!strncmp(tree->GetName(),"TreeTrackReferences",5)) {
232       folder     = "RunMC/Event/Data/TrackReferences";
233       folderName = "TrackReferences" ; 
234     }
235
236     if (!strncmp(tree->GetName(),"TreeD",5)) {
237       folder     = "Run/Event/Data";
238       folderName = "Digits" ; 
239     }
240     if (!strncmp(tree->GetName(),"TreeS",5)) {
241       folder     = "RunMC/Event/Data/SDigits";
242       folderName = "SDigits" ; 
243     }
244     if (!strncmp(tree->GetName(),"TreeR",5)) folder = "Run/Event/RecData";
245
246     if (folder) {
247       if (GetDebug())
248           printf("%15s: Publishing %s to %s\n",ClassName(),name,folder);
249       Publish(folder,address, folderName.Data());
250     }  
251     return branch;
252 }
253
254 //_______________________________________________________________________
255 void AliDetector::Browse(TBrowser *b)
256 {
257   //
258   // Insert Detector objects in the list of objects to be browsed
259   //
260   char name[64];
261   if( fHits == 0) return;
262   TObject *obj;
263   Int_t i, nobjects;
264   //
265   nobjects = fHits->GetEntries();
266   for (i=0;i<nobjects;i++) {
267     obj = fHits->At(i);
268     sprintf(name,"%s_%d",obj->GetName(),i);
269     b->Add(obj, &name[0]);
270   }
271 }
272
273 //_______________________________________________________________________
274 void AliDetector::Copy(AliDetector &) const
275 {
276   //
277   // Copy *this onto det -- not implemented
278   //
279   Fatal("Copy","Not implemented\n");
280 }
281
282 //_______________________________________________________________________
283 void AliDetector::FinishRun()
284 {
285   //
286   // Procedure called at the end of a run.
287   //
288 }
289
290 //_______________________________________________________________________
291 void AliDetector::RemapTrackReferencesIDs(Int_t *map)
292 {
293   // 
294   // Remapping track reference
295   // Called at finish primary
296   //
297   if (!fTrackReferences) return;
298   for (Int_t i=0;i<fTrackReferences->GetEntries();i++){
299     AliTrackReference * ref = dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(i));
300     if (ref) {
301       Int_t newID = map[ref->GetTrack()];
302       if (newID>=0) ref->SetTrack(newID);
303       else ref->SetTrack(-1);
304       
305     }
306   }
307 }
308
309 //_______________________________________________________________________
310 AliHit* AliDetector::FirstHit(Int_t track)
311 {
312   //
313   // Initialise the hit iterator
314   // Return the address of the first hit for track
315   // If track>=0 the track is read from disk
316   // while if track<0 the first hit of the current
317   // track is returned
318   // 
319   if(track>=0) {
320     gAlice->ResetHits();
321     gAlice->TreeH()->GetEvent(track);
322   }
323   //
324   sMaxIterHit=fHits->GetEntriesFast();
325   sCurIterHit=0;
326   if(sMaxIterHit) return dynamic_cast<AliHit*>(fHits->UncheckedAt(0));
327   else            return 0;
328 }
329
330
331 //_______________________________________________________________________
332 AliTrackReference* AliDetector::FirstTrackReference(Int_t track)
333 {
334   //
335   // Initialise the hit iterator
336   // Return the address of the first hit for track
337   // If track>=0 the track is read from disk
338   // while if track<0 the first hit of the current
339   // track is returned
340   // 
341   if(track>=0) {
342     gAlice->ResetTrackReferences();
343     gAlice->TreeTR()->GetEvent(track);
344   }
345   //
346   fMaxIterTrackRef     = fTrackReferences->GetEntriesFast();
347   fCurrentIterTrackRef = 0;
348   if(fMaxIterTrackRef) return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(0));
349   else            return 0;
350 }
351
352 //_______________________________________________________________________
353 AliHit* AliDetector::NextHit()
354 {
355   //
356   // Return the next hit for the current track
357   //
358   if(sMaxIterHit) {
359     if(++sCurIterHit<sMaxIterHit) 
360       return dynamic_cast<AliHit*>(fHits->UncheckedAt(sCurIterHit));
361     else        
362       return 0;
363   } else {
364     printf("* AliDetector::NextHit * Hit Iterator called without calling FistHit before\n");
365     return 0;
366   }
367 }
368
369 //_______________________________________________________________________
370 AliTrackReference* AliDetector::NextTrackReference()
371 {
372   //
373   // Return the next hit for the current track
374   //
375   if(fMaxIterTrackRef) {
376     if(++fCurrentIterTrackRef<fMaxIterTrackRef) 
377       return dynamic_cast<AliTrackReference*>(fTrackReferences->UncheckedAt(fCurrentIterTrackRef));
378     else        
379       return 0;
380   } else {
381     printf("* AliDetector::NextTrackReference * TrackReference  Iterator called without calling FistTrackReference before\n");
382     return 0;
383   }
384 }
385
386 //_______________________________________________________________________
387 void AliDetector::LoadPoints(Int_t)
388 {
389   //
390   // Store x, y, z of all hits in memory
391   //
392   if (fHits == 0) return;
393   //
394   Int_t nhits = fHits->GetEntriesFast();
395   if (nhits == 0) return;
396   Int_t tracks = gAlice->GetNtrack();
397   if (fPoints == 0) fPoints = new TObjArray(tracks);
398   AliHit *ahit;
399   //
400   Int_t *ntrk=new Int_t[tracks];
401   Int_t *limi=new Int_t[tracks];
402   Float_t **coor=new Float_t*[tracks];
403   for(Int_t i=0;i<tracks;i++) {
404     ntrk[i]=0;
405     coor[i]=0;
406     limi[i]=0;
407   }
408   //
409   AliPoints *points = 0;
410   Float_t *fp=0;
411   Int_t trk;
412   Int_t chunk=nhits/4+1;
413   //
414   // Loop over all the hits and store their position
415   for (Int_t hit=0;hit<nhits;hit++) {
416     ahit = dynamic_cast<AliHit*>(fHits->UncheckedAt(hit));
417     trk=ahit->GetTrack();
418     assert(trk<=tracks);
419     if(ntrk[trk]==limi[trk]) {
420       //
421       // Initialise a new track
422       fp=new Float_t[3*(limi[trk]+chunk)];
423       if(coor[trk]) {
424         memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
425         delete [] coor[trk];
426       }
427       limi[trk]+=chunk;
428       coor[trk] = fp;
429     } else {
430       fp = coor[trk];
431     }
432     fp[3*ntrk[trk]  ] = ahit->X();
433     fp[3*ntrk[trk]+1] = ahit->Y();
434     fp[3*ntrk[trk]+2] = ahit->Z();
435     ntrk[trk]++;
436   }
437   //
438   for(trk=0; trk<tracks; ++trk) {
439     if(ntrk[trk]) {
440       points = new AliPoints();
441       points->SetMarkerColor(GetMarkerColor());
442       points->SetMarkerSize(GetMarkerSize());
443       points->SetDetector(this);
444       points->SetParticle(trk);
445       points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
446       fPoints->AddAt(points,trk);
447       delete [] coor[trk];
448       coor[trk]=0;
449     }
450   }
451   delete [] coor;
452   delete [] ntrk;
453   delete [] limi;
454 }
455
456 //_______________________________________________________________________
457 void AliDetector::MakeBranch(Option_t *option, const char *file)
458 {
459   //
460   // Create a new branch in the current Root Tree
461   // The branch of fHits is automatically split
462   //
463  
464   char branchname[10];
465   sprintf(branchname,"%s",GetName());
466   //
467   // Get the pointer to the header
468   const char *cH = strstr(option,"H");
469   //
470   if (fHits && gAlice->TreeH() && cH) {
471     MakeBranchInTree(gAlice->TreeH(), 
472                      branchname, &fHits, fBufferSize, file) ;              
473   }     
474   
475   const char *cD = strstr(option,"D");
476
477   if (cD) {
478     if (file) {
479        fDigitsFile = new char[strlen (file)];
480        strcpy(fDigitsFile,file);
481     }
482   }
483 }
484 //_______________________________________________________________________
485 void AliDetector::MakeBranchTR(Option_t *option, const char *file)
486 {
487   //
488   // Create a new branch in the current Root Tree
489   // The branch of fHits is automatically split
490   //
491  
492   char branchname[10];
493   sprintf(branchname,"%s",GetName());
494   //
495   // Get the pointer to the header
496   const char *cTR = strstr(option,"T");
497   //
498   if (fTrackReferences && gAlice->TreeTR() && cTR) {
499     MakeBranchInTree(gAlice->TreeTR(), 
500                      branchname, &fTrackReferences, fBufferSize, file) ;              
501   }       
502 }
503
504 //_______________________________________________________________________
505 void AliDetector::ResetDigits()
506 {
507   //
508   // Reset number of digits and the digits array
509   //
510   fNdigits   = 0;
511   if (fDigits)   fDigits->Clear();
512 }
513
514 //_______________________________________________________________________
515 void AliDetector::ResetHits()
516 {
517   //
518   // Reset number of hits and the hits array
519   //
520   fNhits   = 0;
521   if (fHits)   fHits->Clear();
522 }
523
524 //_______________________________________________________________________
525 void AliDetector::ResetTrackReferences()
526 {
527   //
528   // Reset number of hits and the hits array
529   //
530   fMaxIterTrackRef   = 0;
531   if (fTrackReferences)   fTrackReferences->Clear();
532 }
533
534 //_______________________________________________________________________
535 void AliDetector::ResetPoints()
536 {
537   //
538   // Reset array of points
539   //
540   if (fPoints) {
541     fPoints->Delete();
542     delete fPoints;
543     fPoints = 0;
544   }
545 }
546
547 //_______________________________________________________________________
548 void AliDetector::SetTreeAddress()
549 {
550   //
551   // Set branch address for the Hits and Digits Trees
552   //
553   TBranch *branch;
554   char branchname[20];
555   sprintf(branchname,"%s",GetName());
556   //
557   // Branch address for hit tree
558   TTree *treeH = gAlice->TreeH();
559   if (treeH && fHits) {
560     branch = treeH->GetBranch(branchname);
561     if (branch) branch->SetAddress(&fHits);
562   }
563   //
564   // Branch address for digit tree
565   TTree *treeD = gAlice->TreeD();
566   if (treeD && fDigits) {
567     branch = treeD->GetBranch(branchname);
568     if (branch) branch->SetAddress(&fDigits);
569   }
570
571   // Branch address for tr  tree
572   TTree *treeTR = gAlice->TreeTR();
573   if (treeTR && fTrackReferences) {
574     branch = treeTR->GetBranch(branchname);
575     if (branch) branch->SetAddress(&fTrackReferences);
576   }
577 }
578
579