]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDetector.cxx
Introducing new Rndm and QA classes
[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 /*
17 $Log$
18 Revision 1.7  2000/10/02 21:28:14  fca
19 Removal of useless dependecies via forward declarations
20
21 Revision 1.6  2000/07/12 08:56:25  fca
22 Coding convention correction and warning removal
23
24 Revision 1.5  1999/09/29 09:24:29  fca
25 Introduction of the Copyright and cvs Log
26
27 */
28
29 ///////////////////////////////////////////////////////////////////////////////
30 //                                                                           //
31 // Base class for ALICE modules. Both sensitive modules (detectors) and      //
32 // non-sensitive ones are described by this base class. This class           //
33 // supports the hit and digit trees produced by the simulation and also      //
34 // the objects produced by the reconstruction.                               //
35 //                                                                           //
36 // This class is also responsible for building the geometry of the           //
37 // detectors.                                                                //
38 //                                                                           //
39 //Begin_Html
40 /*
41 <img src="picts/AliDetectorClass.gif">
42 */
43 //End_Html
44 //                                                                           //
45 ///////////////////////////////////////////////////////////////////////////////
46
47 #include <TTree.h>
48 #include "TBrowser.h"
49
50 #include "AliDetector.h"
51 #include "AliRun.h"
52 #include "AliHit.h"
53 #include "AliPoints.h"
54 // Static variables for the hit iterator routines
55 static Int_t sMaxIterHit=0;
56 static Int_t sCurIterHit=0;
57
58 ClassImp(AliDetector)
59  
60 //_____________________________________________________________________________
61 AliDetector::AliDetector()
62 {
63   //
64   // Default constructor for the AliDetector class
65   //
66   fNhits      = 0;
67   fNdigits    = 0;
68   fPoints     = 0;
69   fHits       = 0;
70   fDigits     = 0;
71   fTimeGate   = 200.e-9;
72   fBufferSize = 16000;
73 }
74  
75 //_____________________________________________________________________________
76 AliDetector::AliDetector(const char* name,const char *title):AliModule(name,title)
77 {
78   //
79   // Normal constructor invoked by all Detectors.
80   // Create the list for detector specific histograms
81   // Add this Detector to the global list of Detectors in Run.
82   //
83
84   fTimeGate   = 200.e-9;
85   fActive     = kTRUE;
86   fNhits      = 0;
87   fHits       = 0;
88   fDigits     = 0;
89   fNdigits    = 0;
90   fPoints     = 0;
91   fBufferSize = 16000;
92 }
93  
94 //_____________________________________________________________________________
95 AliDetector::~AliDetector()
96 {
97   //
98   // Destructor
99   //
100   fNhits      = 0;
101   fNdigits    = 0;
102   //
103   // Delete space point structure
104   if (fPoints) fPoints->Delete();
105   delete fPoints;
106   fPoints     = 0;
107 }
108  
109 //_____________________________________________________________________________
110 void AliDetector::Browse(TBrowser *b)
111 {
112   //
113   // Insert Detector objects in the list of objects to be browsed
114   //
115   char name[64];
116   if( fHits == 0) return;
117   TObject *obj;
118   Int_t i, nobjects;
119   //
120   nobjects = fHits->GetEntries();
121   for (i=0;i<nobjects;i++) {
122     obj = fHits->At(i);
123     sprintf(name,"%s_%d",obj->GetName(),i);
124     b->Add(obj, &name[0]);
125   }
126 }
127
128 //_____________________________________________________________________________
129 void AliDetector::Copy(AliDetector &det) const
130 {
131   //
132   // Copy *this onto det -- not implemented
133   //
134   Fatal("Copy","Not implemented~\n");
135 }
136
137 //_____________________________________________________________________________
138 void AliDetector::FinishRun()
139 {
140   //
141   // Procedure called at the end of a run.
142   //
143 }
144
145 //_____________________________________________________________________________
146 AliHit* AliDetector::FirstHit(Int_t track)
147 {
148   //
149   // Initialise the hit iterator
150   // Return the address of the first hit for track
151   // If track>=0 the track is read from disk
152   // while if track<0 the first hit of the current
153   // track is returned
154   // 
155   if(track>=0) {
156     gAlice->ResetHits();
157     gAlice->TreeH()->GetEvent(track);
158   }
159   //
160   sMaxIterHit=fHits->GetEntriesFast();
161   sCurIterHit=0;
162   if(sMaxIterHit) return (AliHit*) fHits->UncheckedAt(0);
163   else            return 0;
164 }
165
166 //_____________________________________________________________________________
167 AliHit* AliDetector::NextHit()
168 {
169   //
170   // Return the next hit for the current track
171   //
172   if(sMaxIterHit) {
173     if(++sCurIterHit<sMaxIterHit) 
174       return (AliHit*) fHits->UncheckedAt(sCurIterHit);
175     else        
176       return 0;
177   } else {
178     printf("* AliDetector::NextHit * Hit Iterator called without calling FistHit before\n");
179     return 0;
180   }
181 }
182
183 //_____________________________________________________________________________
184 void AliDetector::LoadPoints(Int_t)
185 {
186   //
187   // Store x, y, z of all hits in memory
188   //
189   if (fHits == 0) return;
190   //
191   Int_t nhits = fHits->GetEntriesFast();
192   if (nhits == 0) return;
193   Int_t tracks = gAlice->GetNtrack();
194   if (fPoints == 0) fPoints = new TObjArray(tracks);
195   AliHit *ahit;
196   //
197   Int_t *ntrk=new Int_t[tracks];
198   Int_t *limi=new Int_t[tracks];
199   Float_t **coor=new Float_t*[tracks];
200   for(Int_t i=0;i<tracks;i++) {
201     ntrk[i]=0;
202     coor[i]=0;
203     limi[i]=0;
204   }
205   //
206   AliPoints *points = 0;
207   Float_t *fp=0;
208   Int_t trk;
209   Int_t chunk=nhits/4+1;
210   //
211   // Loop over all the hits and store their position
212   for (Int_t hit=0;hit<nhits;hit++) {
213     ahit = (AliHit*)fHits->UncheckedAt(hit);
214     trk=ahit->GetTrack();
215     if(ntrk[trk]==limi[trk]) {
216       //
217       // Initialise a new track
218       fp=new Float_t[3*(limi[trk]+chunk)];
219       if(coor[trk]) {
220         memcpy(fp,coor[trk],sizeof(Float_t)*3*limi[trk]);
221         delete [] coor[trk];
222       }
223       limi[trk]+=chunk;
224       coor[trk] = fp;
225     } else {
226       fp = coor[trk];
227     }
228     fp[3*ntrk[trk]  ] = ahit->X();
229     fp[3*ntrk[trk]+1] = ahit->Y();
230     fp[3*ntrk[trk]+2] = ahit->Z();
231     ntrk[trk]++;
232   }
233   //
234   for(trk=0; trk<tracks; ++trk) {
235     if(ntrk[trk]) {
236       points = new AliPoints();
237       points->SetMarkerColor(GetMarkerColor());
238       points->SetMarkerSize(GetMarkerSize());
239       points->SetDetector(this);
240       points->SetParticle(trk);
241       points->SetPolyMarker(ntrk[trk],coor[trk],GetMarkerStyle());
242       fPoints->AddAt(points,trk);
243       delete [] coor[trk];
244       coor[trk]=0;
245     }
246   }
247   delete [] coor;
248   delete [] ntrk;
249   delete [] limi;
250 }
251
252 //_____________________________________________________________________________
253 void AliDetector::MakeBranch(Option_t *option)
254 {
255   //
256   // Create a new branch in the current Root Tree
257   // The branch of fHits is automatically split
258   //
259   char branchname[10];
260   sprintf(branchname,"%s",GetName());
261   //
262   // Get the pointer to the header
263   char *cH = strstr(option,"H");
264   //
265   if (fHits   && gAlice->TreeH() && cH) {
266     gAlice->TreeH()->Branch(branchname,&fHits, fBufferSize);
267     printf("* AliDetector::MakeBranch * Making Branch %s for hits\n",branchname);
268   }     
269 }
270
271 //_____________________________________________________________________________
272 void AliDetector::ResetDigits()
273 {
274   //
275   // Reset number of digits and the digits array
276   //
277   fNdigits   = 0;
278   if (fDigits)   fDigits->Clear();
279 }
280
281 //_____________________________________________________________________________
282 void AliDetector::ResetHits()
283 {
284   //
285   // Reset number of hits and the hits array
286   //
287   fNhits   = 0;
288   if (fHits)   fHits->Clear();
289 }
290
291 //_____________________________________________________________________________
292 void AliDetector::ResetPoints()
293 {
294   //
295   // Reset array of points
296   //
297   if (fPoints) {
298     fPoints->Delete();
299     delete fPoints;
300     fPoints = 0;
301   }
302 }
303
304 //_____________________________________________________________________________
305 void AliDetector::SetTreeAddress()
306 {
307   //
308   // Set branch address for the Hits and Digits Trees
309   //
310   TBranch *branch;
311   char branchname[20];
312   sprintf(branchname,"%s",GetName());
313   //
314   // Branch address for hit tree
315   TTree *treeH = gAlice->TreeH();
316   if (treeH && fHits) {
317     branch = treeH->GetBranch(branchname);
318     if (branch) branch->SetAddress(&fHits);
319   }
320   //
321   // Branch address for digit tree
322   TTree *treeD = gAlice->TreeD();
323   if (treeD && fDigits) {
324     branch = treeD->GetBranch(branchname);
325     if (branch) branch->SetAddress(&fDigits);
326   }
327 }
328
329 //_____________________________________________________________________________
330 void AliDetector::Streamer(TBuffer &R__b)
331 {
332   //
333   // Stream an object of class Detector.
334   //
335   if (R__b.IsReading()) {
336     Version_t R__v = R__b.ReadVersion(); if (R__v) { }
337     TNamed::Streamer(R__b);
338     TAttLine::Streamer(R__b);
339     TAttMarker::Streamer(R__b);
340     AliModule::Streamer(R__b);
341     R__b >> fTimeGate;
342     R__b >> fIshunt;
343     //R__b >> fNhits;
344     //
345     // Stream the pointers but not the TClonesArrays
346     R__b >> fHits; // diff
347     R__b >> fDigits; // diff
348     
349   } else {
350     R__b.WriteVersion(AliDetector::IsA());
351     TNamed::Streamer(R__b);
352     TAttLine::Streamer(R__b);
353     TAttMarker::Streamer(R__b);
354     AliModule::Streamer(R__b);
355     R__b << fTimeGate;
356     R__b << fIshunt;
357     //R__b << fNhits;
358     //
359     // Stream the pointers but not the TClonesArrays
360     R__b << fHits; // diff
361     R__b << fDigits; // diff
362   }
363 }
364