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