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