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