]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliDetector.cxx
Do not save CVS subdirectories
[u/mrichter/AliRoot.git] / STEER / AliDetector.cxx
CommitLineData
fe4da5cc 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 "AliMC.h"
23#include <TClass.h>
24#include <TNode.h>
25#include <TRandom.h>
26
27// Static variables for the hit iterator routines
28static Int_t sMaxIterHit=0;
29static Int_t sCurIterHit=0;
30
31ClassImp(AliDetector)
32
33//_____________________________________________________________________________
34AliDetector::AliDetector()
35{
36 //
37 // Default constructor for the AliDetector class
38 //
39 fNhits = 0;
40 fNdigits = 0;
41 fHistograms = 0;
42 fNodes = 0;
43 fPoints = 0;
44 fHits = 0;
45 fDigits = 0;
46 fTimeGate = 200.e-9;
47 fActive = kTRUE;
48 fBufferSize = 16000;
49}
50
51//_____________________________________________________________________________
52AliDetector::AliDetector(const char* name,const char *title):TNamed(name,title)
53{
54 //
55 // Normal constructor invoked by all Detectors.
56 // Create the list for detector specific histograms
57 // Add this Detector to the global list of Detectors in Run.
58 //
59
60 fTimeGate = 200.e-9;
61 fActive = kTRUE;
62 fNhits = 0;
63 fHits = 0;
64 fDigits = 0;
65 fNdigits = 0;
66 fPoints = 0;
67 fBufferSize = 16000;
68 //
69 // Initialises the histogram list
70 fHistograms = new TList();
71 //
72 // Initialises the list of ROOT TNodes
73 fNodes = new TList();
74 //
75 // Get the detector numeric ID
76 Int_t id = gAlice->GetDetectorID(name);
77 if (id < 0) {
78 // Unknown detector !
79 printf(" * AliRun::Ctor * ERROR Unknown detector: %s\n",name);
80 return;
81 }
82 //
83 // Add this detector to the list of detectors
84 gAlice->Detectors()->AddAtAndExpand(this,id);
85 //
86 //
87 SetMarkerColor(3);
88 //
89 // Allocate space for tracking media and material indexes
90 fIdtmed = new TArrayI(100);
91 fIdmate = new TArrayI(100);
92 for(Int_t i=0;i<100;i++) (*fIdmate)[i]=(*fIdtmed)[i]=0;
93 //
94 // Prepare to find the tracking media range
95 fLoMedium = 65536;
96 fHiMedium = 0;
97}
98
99//_____________________________________________________________________________
100AliDetector::~AliDetector()
101{
102 //
103 // Destructor
104 //
105 fNhits = 0;
106 fNdigits = 0;
107 fHistograms = 0;
108 //
109 // Delete ROOT geometry
110 fNodes->Clear();
111 delete fNodes;
112 //
113 // Delete space point structure
114 if (fPoints) fPoints->Delete();
115 delete fPoints;
116 fPoints = 0;
117 //
118 // Delete TArray objects
119 delete fIdtmed;
120 delete fIdmate;
121}
122
123//_____________________________________________________________________________
124void AliDetector::Browse(TBrowser *b)
125{
126 //
127 // Insert Detector objects in the list of objects to be browsed
128 //
129 char name[64];
130 if( fHits == 0) return;
131 TObject *obj;
132 Int_t i, nobjects;
133 //
134 nobjects = fHits->GetEntries();
135 for (i=0;i<nobjects;i++) {
136 obj = fHits->At(i);
137 sprintf(name,"%s_%d",obj->GetName(),i);
138 b->Add(obj, &name[0]);
139 }
140}
141
142//_____________________________________________________________________________
143void AliDetector::Disable()
144{
145 //
146 // Disable detector on viewer
147 //
148 fActive = kFALSE;
149 TIter next(fNodes);
150 TNode *node;
151 //
152 // Loop through geometry to disable all
153 // nodes for this detector
154 while((node = (TNode*)next())) {
155 node->SetVisibility(0);
156 }
157}
158
159//_____________________________________________________________________________
160Int_t AliDetector::DistancetoPrimitive(Int_t, Int_t)
161{
162 //
163 // Return distance from mouse pointer to object
164 // Dummy routine for the moment
165 //
166 return 9999;
167}
168
169//_____________________________________________________________________________
170void AliDetector::Enable()
171{
172 //
173 // Enable detector on the viewver
174 //
175 fActive = kTRUE;
176 TIter next(fNodes);
177 TNode *node;
178 //
179 // Loop through geometry to enable all
180 // nodes for this detector
181 while((node = (TNode*)next())) {
182 node->SetVisibility(1);
183 }
184}
185
186//_____________________________________________________________________________
187void AliDetector::FinishRun()
188{
189 //
190 // Procedure called at the end of a run.
191 //
192}
193
194//_____________________________________________________________________________
195AliHit* AliDetector::FirstHit(Int_t track)
196{
197 //
198 // Initialise the hit iterator
199 // Return the address of the first hit for track
200 // If track>=0 the track is read from disk
201 // while if track<0 the first hit of the current
202 // track is returned
203 //
204 if(track>=0) {
205 gAlice->ResetHits();
206 gAlice->TreeH()->GetEvent(track);
207 }
208 //
209 sMaxIterHit=fHits->GetEntriesFast();
210 sCurIterHit=0;
211 if(sMaxIterHit) return (AliHit*) fHits->UncheckedAt(0);
212 else return 0;
213}
214
215//_____________________________________________________________________________
216AliHit* AliDetector::NextHit()
217{
218 //
219 // Return the next hit for the current track
220 //
221 if(sMaxIterHit) {
222 if(++sCurIterHit<sMaxIterHit)
223 return (AliHit*) fHits->UncheckedAt(sCurIterHit);
224 else
225 return 0;
226 } else {
227 printf("* AliDetector::NextHit * Hit Iterator called without calling FistHit before\n");
228 return 0;
229 }
230}
231
232//_____________________________________________________________________________
233void AliDetector::LoadPoints(Int_t)
234{
235 //
236 // Store x, y, z of all hits in memory
237 //
238 if (fHits == 0) return;
239 //
240 if (fPoints == 0) fPoints = new TObjArray(gAlice->GetNtrack());
241 Int_t nhits = fHits->GetEntriesFast();
242 if (nhits == 0) return;
243 AliHit *ahit;
244 //
245 AliPoints *points = 0;
246 Int_t trko=-99, trk;
247 //
248 // Loop over all the hits and store their position
249 for (Int_t hit=0;hit<nhits;hit++) {
250 ahit = (AliHit*)fHits->UncheckedAt(hit);
251 if(trko!=(trk=ahit->GetTrack())) {
252 //
253 // Initialise a new track
254 trko=trk;
255 points = new AliPoints(nhits);
256 fPoints->AddAt(points,trk);
257 points->SetMarkerColor(GetMarkerColor());
258 points->SetMarkerStyle(GetMarkerStyle());
259 points->SetMarkerSize(GetMarkerSize());
260 points->SetDetector(this);
261 points->SetParticle(trk);
262 }
263 points->SetPoint(hit,ahit->fX,ahit->fY,ahit->fZ);
264 }
265}
266
267//_____________________________________________________________________________
268void AliDetector::MakeBranch(Option_t *option)
269{
270 //
271 // Create a new branch in the current Root Tree
272 // The branch of fHits is automatically split
273 //
274 char branchname[10];
275 sprintf(branchname,"%s",GetName());
276 //
277 // Get the pointer to the header
278 char *H = strstr(option,"H");
279 //
280 if (fHits && gAlice->TreeH() && H) {
281 gAlice->TreeH()->Branch(branchname,&fHits, fBufferSize);
282 printf("* AliDetector::MakeBranch * Making Branch %s for hits\n",branchname);
283 }
284}
285
286//_____________________________________________________________________________
287void AliDetector::AliMaterial(Int_t imat, const char* name, Float_t a,
288 Float_t z, Float_t dens, Float_t radl,
289 Float_t absl, Float_t *buf, Int_t nwbuf) const
290{
291 //
292 // Store the parameters for a material
293 //
294 // imat the material index will be stored in (*fIdmate)[imat]
295 // name material name
296 // a atomic mass
297 // z atomic number
298 // dens density
299 // radl radiation length
300 // absl absorbtion length
301 // buf adress of an array user words
302 // nwbuf number of user words
303 //
304 Int_t kmat;
305 AliMC::GetMC()->Material(kmat, name, a, z, dens, radl, absl, buf, nwbuf);
306 (*fIdmate)[imat]=kmat;
307}
308
309
310//_____________________________________________________________________________
311void AliDetector::AliMixture(Int_t imat, const char *name, Float_t *a,
312 Float_t *z, Float_t dens, Int_t nlmat,
313 Float_t *wmat) const
314{
315 //
316 // Defines mixture or compound imat as composed by
317 // nlmat materials defined by arrays a, z and wmat
318 //
319 // If nlmat > 0 wmat contains the proportion by
320 // weights of each basic material in the mixture
321 //
322 // If nlmat < 0 wmat contains the number of atoms
323 // of eack kind in the molecule of the compound
324 // In this case, wmat is changed on output to the relative weigths.
325 //
326 // imat the material index will be stored in (*fIdmate)[imat]
327 // name material name
328 // a array of atomic masses
329 // z array of atomic numbers
330 // dens density
331 // nlmat number of components
332 // wmat array of concentrations
333 //
334 Int_t kmat;
335 AliMC::GetMC()->Mixture(kmat, name, a, z, dens, nlmat, wmat);
336 (*fIdmate)[imat]=kmat;
337}
338
339//_____________________________________________________________________________
340void AliDetector::AliMedium(Int_t numed, const char *name, Int_t nmat,
341 Int_t isvol, Int_t ifield, Float_t fieldm,
342 Float_t tmaxfd, Float_t stemax, Float_t deemax,
343 Float_t epsil, Float_t stmin, Float_t *ubuf,
344 Int_t nbuf) const
345{
346 //
347 // Store the parameters of a tracking medium
348 //
349 // numed the medium number is stored into (*fIdtmed)[numed-1]
350 // name medium name
351 // nmat the material number is stored into (*fIdmate)[nmat]
352 // isvol sensitive volume if isvol!=0
353 // ifield magnetic field flag (see below)
354 // fieldm maximum magnetic field
355 // tmaxfd maximum deflection angle due to magnetic field
356 // stemax maximum step allowed
357 // deemax maximum fractional energy loss in one step
358 // epsil tracking precision in cm
359 // stmin minimum step due to continuous processes
360 //
361 // ifield = 0 no magnetic field
362 // = -1 user decision in guswim
363 // = 1 tracking performed with Runge Kutta
364 // = 2 tracking performed with helix
365 // = 3 constant magnetic field along z
366 //
367 Int_t kmed;
368 Int_t *idtmed = gAlice->Idtmed();
369 AliMC::GetMC()->Medium(kmed,name, (*fIdmate)[nmat], isvol, ifield, fieldm,
370 tmaxfd, stemax, deemax, epsil, stmin, ubuf, nbuf);
371 idtmed[numed-1]=kmed;
372}
373
374//_____________________________________________________________________________
375void AliDetector::AliMatrix(Int_t &nmat, Float_t theta1, Float_t phi1,
376 Float_t theta2, Float_t phi2, Float_t theta3,
377 Float_t phi3) const
378{
379 //
380 // Define a rotation matrix. Angles are in degrees.
381 //
382 // nmat on output contains the number assigned to the rotation matrix
383 // theta1 polar angle for axis I
384 // phi1 azimuthal angle for axis I
385 // theta2 polar angle for axis II
386 // phi2 azimuthal angle for axis II
387 // theta3 polar angle for axis III
388 // phi3 azimuthal angle for axis III
389 //
390 AliMC::GetMC()->Matrix(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
391}
392
393//_____________________________________________________________________________
394void AliDetector::ResetDigits()
395{
396 //
397 // Reset number of digits and the digits array
398 //
399 fNdigits = 0;
400 if (fDigits) fDigits->Clear();
401}
402
403//_____________________________________________________________________________
404void AliDetector::ResetHits()
405{
406 //
407 // Reset number of hits and the hits array
408 //
409 fNhits = 0;
410 if (fHits) fHits->Clear();
411}
412
413//_____________________________________________________________________________
414void AliDetector::ResetPoints()
415{
416 //
417 // Reset array of points
418 //
419 if (fPoints) {
420 fPoints->Delete();
421 delete fPoints;
422 fPoints = 0;
423 }
424}
425
426
427//_____________________________________________________________________________
428void AliDetector::StepManager()
429{
430 //
431 // Procedure called at every step inside the detector
432 //
433 printf("* AliDetector::StepManager * Generic Step Manager called for Detector: %s\n",fName.Data());
434}
435
436//_____________________________________________________________________________
437void AliDetector::SetEuclidFile(char* material, char* geometry)
438{
439 //
440 // Sets the name of the Euclid file
441 //
442 fEuclidMaterial=material;
443 if(geometry) {
444 fEuclidGeometry=geometry;
445 } else {
446 char* name = new char[strlen(material)];
447 strcpy(name,material);
448 strcpy(&name[strlen(name)-4],".euc");
449 fEuclidGeometry=name;
450 delete [] name;
451 }
452}
453
454//_____________________________________________________________________________
455void AliDetector::SetTreeAddress()
456{
457 //
458 // Set branch address for the Hits and Digits Trees
459 //
460 TBranch *branch;
461 char branchname[20];
462 sprintf(branchname,"%s",GetName());
463 //
464 // Branch address for hit tree
465 TTree *treeH = gAlice->TreeH();
466 if (treeH && fHits) {
467 branch = treeH->GetBranch(branchname);
468 if (branch) branch->SetAddress(&fHits);
469 }
470 //
471 // Branch address for digit tree
472 TTree *treeD = gAlice->TreeD();
473 if (treeD && fDigits) {
474 branch = treeD->GetBranch(branchname);
475 if (branch) branch->SetAddress(&fDigits);
476 }
477}
478
479//_____________________________________________________________________________
480void AliDetector::Streamer(TBuffer &R__b)
481{
482 //
483 // Stream an object of class Detector.
484 //
485 if (R__b.IsReading()) {
486 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
487 TNamed::Streamer(R__b);
488 TAttLine::Streamer(R__b);
489 TAttMarker::Streamer(R__b);
490 fEuclidMaterial.Streamer(R__b);
491 fEuclidGeometry.Streamer(R__b);
492 R__b >> fTimeGate;
493 R__b >> fActive;
494 R__b >> fIshunt;
495 //R__b >> fNhits;
496 R__b >> fHistograms;
497 //
498 // Stream the pointers but not the TClonesArrays
499 R__b >> fNodes; // diff
500
501 R__b >> fHits; // diff
502 R__b >> fDigits; // diff
503
504 } else {
505 R__b.WriteVersion(AliDetector::IsA());
506 TNamed::Streamer(R__b);
507 TAttLine::Streamer(R__b);
508 TAttMarker::Streamer(R__b);
509 fEuclidMaterial.Streamer(R__b);
510 fEuclidGeometry.Streamer(R__b);
511 R__b << fTimeGate;
512 R__b << fActive;
513 R__b << fIshunt;
514 //R__b << fNhits;
515 R__b << fHistograms;
516 //
517 // Stream the pointers but not the TClonesArrays
518 R__b << fNodes; // diff
519
520 R__b << fHits; // diff
521 R__b << fDigits; // diff
522 }
523}
524