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