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