]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliDetector.cxx
Put headers before libraries in the make.
[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="gif/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   if (fPoints == 0) fPoints = new TObjArray(gAlice->GetNtrack());
155   Int_t nhits = fHits->GetEntriesFast();
156   if (nhits == 0) return;
157   AliHit *ahit;
158   //
159   AliPoints *points = 0;
160   Int_t trko=-99, trk;
161   //
162   // Loop over all the hits and store their position
163   for (Int_t hit=0;hit<nhits;hit++) {
164     ahit = (AliHit*)fHits->UncheckedAt(hit);
165     if(trko!=(trk=ahit->GetTrack())) {
166       //
167       // Initialise a new track
168       trko=trk;
169       points = new AliPoints(nhits);
170       fPoints->AddAt(points,trk);
171       points->SetMarkerColor(GetMarkerColor());
172       points->SetMarkerStyle(GetMarkerStyle());
173       points->SetMarkerSize(GetMarkerSize());
174       points->SetDetector(this);
175       points->SetParticle(trk);
176     }
177     points->SetPoint(hit,ahit->fX,ahit->fY,ahit->fZ);
178   }
179 }
180
181 //_____________________________________________________________________________
182 void AliDetector::MakeBranch(Option_t *option)
183 {
184   //
185   // Create a new branch in the current Root Tree
186   // The branch of fHits is automatically split
187   //
188   char branchname[10];
189   sprintf(branchname,"%s",GetName());
190   //
191   // Get the pointer to the header
192   char *H = strstr(option,"H");
193   //
194   if (fHits   && gAlice->TreeH() && H) {
195     gAlice->TreeH()->Branch(branchname,&fHits, fBufferSize);
196     printf("* AliDetector::MakeBranch * Making Branch %s for hits\n",branchname);
197   }     
198 }
199
200 //_____________________________________________________________________________
201 void AliDetector::ResetDigits()
202 {
203   //
204   // Reset number of digits and the digits array
205   //
206   fNdigits   = 0;
207   if (fDigits)   fDigits->Clear();
208 }
209
210 //_____________________________________________________________________________
211 void AliDetector::ResetHits()
212 {
213   //
214   // Reset number of hits and the hits array
215   //
216   fNhits   = 0;
217   if (fHits)   fHits->Clear();
218 }
219
220 //_____________________________________________________________________________
221 void AliDetector::ResetPoints()
222 {
223   //
224   // Reset array of points
225   //
226   if (fPoints) {
227     fPoints->Delete();
228     delete fPoints;
229     fPoints = 0;
230   }
231 }
232
233 //_____________________________________________________________________________
234 void AliDetector::SetTreeAddress()
235 {
236   //
237   // Set branch address for the Hits and Digits Trees
238   //
239   TBranch *branch;
240   char branchname[20];
241   sprintf(branchname,"%s",GetName());
242   //
243   // Branch address for hit tree
244   TTree *treeH = gAlice->TreeH();
245   if (treeH && fHits) {
246     branch = treeH->GetBranch(branchname);
247     if (branch) branch->SetAddress(&fHits);
248   }
249   //
250   // Branch address for digit tree
251   TTree *treeD = gAlice->TreeD();
252   if (treeD && fDigits) {
253     branch = treeD->GetBranch(branchname);
254     if (branch) branch->SetAddress(&fDigits);
255   }
256 }
257
258 //_____________________________________________________________________________
259 void AliDetector::Streamer(TBuffer &R__b)
260 {
261   //
262   // Stream an object of class Detector.
263   //
264   if (R__b.IsReading()) {
265     Version_t R__v = R__b.ReadVersion(); if (R__v) { }
266     TNamed::Streamer(R__b);
267     TAttLine::Streamer(R__b);
268     TAttMarker::Streamer(R__b);
269     AliModule::Streamer(R__b);
270     R__b >> fTimeGate;
271     R__b >> fIshunt;
272     //R__b >> fNhits;
273     //
274     // Stream the pointers but not the TClonesArrays
275     R__b >> fHits; // diff
276     R__b >> fDigits; // diff
277     
278   } else {
279     R__b.WriteVersion(AliDetector::IsA());
280     TNamed::Streamer(R__b);
281     TAttLine::Streamer(R__b);
282     TAttMarker::Streamer(R__b);
283     AliModule::Streamer(R__b);
284     R__b << fTimeGate;
285     R__b << fIshunt;
286     //R__b << fNhits;
287     //
288     // Stream the pointers but not the TClonesArrays
289     R__b << fHits; // diff
290     R__b << fDigits; // diff
291   }
292 }
293