]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRD.cxx
Geometry update, Removal of compiler warnings
[u/mrichter/AliRoot.git] / TRD / AliTRD.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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 //  Transition Radiation Detector                                            //
21 //  This class contains the basic functions for the Transition Radiation     //
22 //  Detector.                                                                //
23 //                                                                           //
24 ///////////////////////////////////////////////////////////////////////////////
25
26 #include <stdlib.h>
27 #include <Riostream.h>
28
29 #include <TFile.h>
30 #include <TGeometry.h>
31 #include <TLorentzVector.h>
32 #include <TMath.h>
33 #include <TNode.h>
34 #include <TPGON.h> 
35 #include <TParticle.h>
36 #include <TROOT.h>
37 #include <TTree.h>
38 #include <TVirtualMC.h>
39  
40 #include "AliConst.h"
41 #include "AliDigit.h"
42 #include "AliLoader.h"
43 #include "AliMagF.h"
44 #include "AliRun.h"
45 #include "AliTRD.h"
46 #include "AliTRDcluster.h"
47 #include "AliTRDclusterizer.h"
48 #include "AliTRDdigit.h"
49 #include "AliTRDdigitizer.h"
50 #include "AliTRDdigitsManager.h"
51 #include "AliTRDgeometryFull.h"
52 #include "AliTRDgeometryHole.h"
53 #include "AliTRDhit.h"
54 #include "AliTRDpoints.h"
55 #include "AliTRDrecPoint.h"
56 #include "AliTRDtrackHits.h"  
57 #include "AliTrackReference.h"
58
59 ClassImp(AliTRD)
60  
61 //_____________________________________________________________________________
62 AliTRD::AliTRD()
63 {
64   //
65   // Default constructor
66   //
67
68   fIshunt        = 0;
69   fGasMix        = 0;
70   fHits          = 0;
71   fDigits        = 0;
72
73   fRecPoints     = 0;
74   fNRecPoints    = 0;
75
76   fGeometry      = 0;
77
78   fGasDensity    = 0;
79   fFoilDensity   = 0;
80
81   fDrawTR        = 0;
82   fDisplayType   = 0;
83  
84   fTrackHits     = 0; 
85   fHitType       = 0; 
86
87 }
88  
89 //_____________________________________________________________________________
90 AliTRD::AliTRD(const char *name, const char *title)
91        : AliDetector(name,title)
92 {
93   //
94   // Standard constructor for the TRD
95   //
96
97   // Check that FRAME is there otherwise we have no place where to
98   // put TRD
99   AliModule* frame = gAlice->GetModule("FRAME");
100   if (!frame) {
101     Error("Ctor","TRD needs FRAME to be present\n");
102     exit(1);
103   } 
104
105   // Define the TRD geometry according to the FRAME geometry
106   if      (frame->IsVersion() == 0) {
107     // Geometry with hole
108     fGeometry = new AliTRDgeometryHole();
109   }
110   else if (frame->IsVersion() == 1) {
111     // Geometry without hole
112     fGeometry = new AliTRDgeometryFull();
113   }
114   else {
115     Error("Ctor","Could not find valid FRAME version\n");
116     exit(1);
117   }
118
119   // Allocate the hit array
120   fHits           = new TClonesArray("AliTRDhit"     ,405);
121   gAlice->AddHitList(fHits);
122
123   // Allocate the digits array
124   fDigits         = 0;
125
126   // Allocate the rec point array
127   fRecPoints     = new TObjArray(400);
128   fNRecPoints    = 0;
129    
130   fIshunt        = 0;
131   fGasMix        = 1;
132
133   fGasDensity    = 0;
134   fFoilDensity   = 0;
135
136   fDrawTR        = 0;
137   fDisplayType   = 0;
138
139   fTrackHits     = 0;
140   fHitType       = 2;
141
142   SetMarkerColor(kWhite);   
143
144 }
145
146 //_____________________________________________________________________________
147 AliTRD::AliTRD(const AliTRD &trd):AliDetector(trd)
148 {
149   //
150   // Copy constructor
151   //
152
153   ((AliTRD &) trd).Copy(*this);
154
155 }
156
157 //_____________________________________________________________________________
158 AliTRD::~AliTRD()
159 {
160   //
161   // TRD destructor
162   //
163
164   fIshunt = 0;
165
166   if (fGeometry) {
167     delete fGeometry;
168     fGeometry  = 0;
169   }
170   if (fHits) {
171     delete fHits;
172     fHits      = 0;
173   }
174   if (fRecPoints) {
175     delete fRecPoints;
176     fRecPoints = 0;
177   }
178   if (fTrackHits) {
179     delete fTrackHits;
180     fTrackHits = 0;
181   }
182
183 }
184
185 //_____________________________________________________________________________
186 void AliTRD::AddCluster(Float_t *pos, Int_t det, Float_t amp
187                       , Int_t *tracks, Float_t *sig, Int_t iType)
188 {
189   //
190   // Add a cluster for the TRD
191   //
192
193   AliTRDcluster *c = new AliTRDcluster();
194
195   c->SetDetector(det);
196   c->AddTrackIndex(tracks);
197   c->SetQ(amp);
198   c->SetY(pos[0]);
199   c->SetZ(pos[1]);
200   c->SetSigmaY2(sig[0]);   
201   c->SetSigmaZ2(sig[1]);
202   c->SetLocalTimeBin(((Int_t) pos[2]));
203
204   switch (iType) {
205   case 0:
206     c->Set2pad();
207     break;
208   case 1:
209     c->Set3pad();
210     break;
211   case 2:
212     c->Set4pad();
213     break;
214   case 3:
215     c->Set5pad();
216     break;
217   case 4:
218     c->SetLarge();
219     break;
220   };
221
222   fRecPoints->Add(c);
223
224 }
225
226 //_____________________________________________________________________________
227 void AliTRD::Hits2Digits()
228 {
229   //
230   // Create digits
231   //
232   if (!fLoader->TreeH()) fLoader->LoadHits("read");
233   SetTreeAddress();
234   
235   AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
236                                                   ,"TRD digitizer class");
237   digitizer->SetDebug(GetDebug());
238   digitizer->SetEvent(gAlice->GetEvNumber());
239   
240   // Initialization
241   digitizer->InitDetector();
242     
243   // Create the digits
244   digitizer->MakeDigits();
245   
246   // Write the digits into the input file
247
248   if (!fLoader->TreeD()) fLoader->MakeTree("D");
249   if (digitizer->MakeBranch(fLoader->TreeD())) {
250
251     digitizer->WriteDigits();
252
253     // Save the digitizer class in the AliROOT 
254     digitizer->Write();
255
256   }
257
258 }
259
260 //_____________________________________________________________________________
261 void AliTRD::Hits2SDigits()
262 {
263   //
264   // Create summable digits
265   //
266   if (!fLoader->TreeH()) fLoader->LoadHits("read");
267   SetTreeAddress();
268
269   AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
270                                                   ,"TRD digitizer class");
271   digitizer->SetDebug(GetDebug());
272
273   // For the summable digits
274   digitizer->SetSDigits(kTRUE);
275   digitizer->SetEvent(gAlice->GetEvNumber());
276
277   // Initialization
278   digitizer->InitDetector();
279     
280   // Create the TRD s-digits branch
281   digitizer->MakeDigits();
282   
283   // Write the digits into the input file
284   if (!fLoader->TreeS()) fLoader->MakeTree("S");
285   
286   if (digitizer->MakeBranch(fLoader->TreeS())) {
287
288     digitizer->WriteDigits();
289
290     // Save the digitizer class in the AliROOT 
291     digitizer->Write();
292
293   }
294
295 }
296
297 //_____________________________________________________________________________
298 void AliTRD::SDigits2Digits()
299 {
300   //
301   // Create final digits from summable digits
302   //
303
304    // Create the TRD digitizer
305   AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
306                                                   ,"TRD digitizer class");  
307   digitizer->SetDebug(GetDebug());
308
309   // Set the parameter
310   digitizer->SetEvent(gAlice->GetEvNumber());
311
312   // Initialization
313   digitizer->InitDetector();
314
315   // Read the s-digits via digits manager
316   AliTRDdigitsManager *sdigitsManager = new AliTRDdigitsManager();
317  
318   sdigitsManager->SetDebug(GetDebug());
319   sdigitsManager->SetSDigits(kTRUE);
320   sdigitsManager->CreateArrays();
321   
322   if (!fLoader->TreeS()) 
323     if (fLoader->LoadSDigits("read"))
324      {
325        Error("SDigits2Digits","Error while reading SDigits for event %d",gAlice->GetEvNumber());
326        return;
327      }
328   if (!fLoader->TreeS()) return;
329   
330   sdigitsManager->ReadDigits(fLoader->TreeS());
331
332   // Add the s-digits to the input list 
333   digitizer->AddSDigitsManager(sdigitsManager);
334
335   // Convert the s-digits to normal digits
336   digitizer->SDigits2Digits();
337
338   // Store the digits
339   if (!fLoader->TreeD()) fLoader->MakeTree("D");
340   if (digitizer->MakeBranch(fLoader->TreeD())){
341     digitizer->WriteDigits();
342   }
343
344 }
345
346 //_____________________________________________________________________________
347 void AliTRD::AddHit(Int_t track, Int_t det, Float_t *hits, Int_t q
348                   , Bool_t inDrift)
349 {
350   //
351   // Add a hit for the TRD
352   // 
353   // The data structure is set according to fHitType:
354   //   bit0: standard TClonesArray
355   //   bit1: compressed trackHits structure
356   //
357
358   if (fHitType & 1) {
359     TClonesArray &lhits = *fHits;
360     new(lhits[fNhits++]) AliTRDhit(fIshunt,track,det,hits,q);
361   }
362
363   if (fHitType > 1) {
364     AddHit2(track,det,hits,q,inDrift);
365   }
366
367 }
368
369 //_____________________________________________________________________________
370 void AliTRD::BuildGeometry()
371 {
372   //
373   // Create the ROOT TNode geometry for the TRD
374   //
375
376   TNode *node, *top;
377   TPGON *pgon;
378
379   Float_t rmin, rmax;
380   Float_t zmax1, zmax2;
381
382   Int_t   iPlan;
383  
384   const Int_t kColorTRD = 46;
385   
386   // Find the top node alice
387   top = gAlice->GetGeometry()->GetNode("alice");
388   
389   if      (fDisplayType == 0) {
390
391     pgon = new TPGON("S_TRD","TRD","void",0,360,AliTRDgeometry::Nsect(),4);
392     rmin = AliTRDgeometry::Rmin();
393     rmax = AliTRDgeometry::Rmax();
394     pgon->DefineSection(0,-AliTRDgeometry::Zmax1(),rmax,rmax);
395     pgon->DefineSection(1,-AliTRDgeometry::Zmax2(),rmin,rmax);
396     pgon->DefineSection(2, AliTRDgeometry::Zmax2(),rmin,rmax);
397     pgon->DefineSection(3, AliTRDgeometry::Zmax1(),rmax,rmax);
398     top->cd();
399     node = new TNode("TRD","TRD","S_TRD",0,0,0,"");
400     node->SetLineColor(kColorTRD);
401     fNodes->Add(node);
402
403   }
404   else if (fDisplayType == 1) {
405
406     Char_t name[7];
407
408     Float_t slope = (AliTRDgeometry::Zmax1() - AliTRDgeometry::Zmax2())
409                   / (AliTRDgeometry::Rmax()  - AliTRDgeometry::Rmin());
410
411     rmin  = AliTRDgeometry::Rmin() + AliTRDgeometry::CraHght();
412     rmax  = rmin                   + AliTRDgeometry::CdrHght();
413
414     Float_t thickness = rmin - AliTRDgeometry::Rmin();
415     zmax2 = AliTRDgeometry::Zmax2() + slope * thickness;
416     zmax1 = zmax2 + slope * AliTRDgeometry::DrThick();
417
418     for (iPlan = 0; iPlan < AliTRDgeometry::Nplan(); iPlan++) {
419
420       sprintf(name,"S_TR1%d",iPlan);
421       pgon  = new TPGON(name,"TRD","void",0,360,AliTRDgeometry::Nsect(),4);
422       pgon->DefineSection(0,-zmax1,rmax,rmax);
423       pgon->DefineSection(1,-zmax2,rmin,rmax);
424       pgon->DefineSection(2, zmax2,rmin,rmax);
425       pgon->DefineSection(3, zmax1,rmax,rmax);
426       top->cd();
427       node = new TNode("TRD","TRD",name,0,0,0,"");
428       node->SetLineColor(kColorTRD);
429       fNodes->Add(node);
430
431       Float_t height = AliTRDgeometry::Cheight() + AliTRDgeometry::Cspace(); 
432       rmin  = rmin  + height;
433       rmax  = rmax  + height;
434       zmax1 = zmax1 + slope * height;
435       zmax2 = zmax2 + slope * height;
436
437     }
438
439     thickness += AliTRDgeometry::DrThick();
440     rmin  = AliTRDgeometry::Rmin() + thickness;
441     rmax  = rmin + AliTRDgeometry::AmThick();
442     zmax2 = AliTRDgeometry::Zmax2() + slope * thickness;
443     zmax1 = zmax2 + slope * AliTRDgeometry::AmThick();
444
445     for (iPlan = 0; iPlan < AliTRDgeometry::Nplan(); iPlan++) {
446
447       sprintf(name,"S_TR2%d",iPlan);
448       pgon  = new TPGON(name,"TRD","void",0,360,AliTRDgeometry::Nsect(),4);
449       pgon->DefineSection(0,-zmax1,rmax,rmax);
450       pgon->DefineSection(1,-zmax2,rmin,rmax);
451       pgon->DefineSection(2, zmax2,rmin,rmax);
452       pgon->DefineSection(3, zmax1,rmax,rmax);
453       top->cd();
454       node = new TNode("TRD","TRD",name,0,0,0,"");
455       node->SetLineColor(kColorTRD);
456       fNodes->Add(node);
457
458       Float_t height = AliTRDgeometry::Cheight() + AliTRDgeometry::Cspace(); 
459       rmin  = rmin  + height;
460       rmax  = rmax  + height;
461       zmax1 = zmax1 + slope * height;
462       zmax2 = zmax2 + slope * height;
463
464     }
465
466   }
467
468 }
469  
470 //_____________________________________________________________________________
471 void AliTRD::Copy(TObject &trd)
472 {
473   //
474   // Copy function
475   //
476
477   ((AliTRD &) trd).fGasMix      = fGasMix;
478   ((AliTRD &) trd).fGeometry    = fGeometry;       
479   ((AliTRD &) trd).fRecPoints   = fRecPoints;
480   ((AliTRD &) trd).fNRecPoints  = fNRecPoints;
481   ((AliTRD &) trd).fGasDensity  = fGasDensity;
482   ((AliTRD &) trd).fFoilDensity = fFoilDensity;
483   ((AliTRD &) trd).fDrawTR      = fDrawTR;
484   ((AliTRD &) trd).fDisplayType = fDisplayType;
485   ((AliTRD &) trd).fHitType     = fHitType;
486
487   //AliDetector::Copy(trd);
488
489 }
490
491 //_____________________________________________________________________________
492 void AliTRD::CreateGeometry()
493 {
494   //
495   // Creates the volumes for the TRD chambers
496   //
497
498   // Check that FRAME is there otherwise we have no place where to put the TRD
499   AliModule* frame = gAlice->GetModule("FRAME");
500   if (!frame) {
501     printf(" The TRD needs the FRAME to be defined first\n");
502     return;
503   }
504
505   fGeometry->CreateGeometry(fIdtmed->GetArray() - 1299);
506
507 }
508  
509 //_____________________________________________________________________________
510 void AliTRD::CreateMaterials()
511 {
512   //
513   // Create the materials for the TRD
514   // Origin Y.Foka
515   //
516
517   Int_t   isxfld = gAlice->Field()->Integ();
518   Float_t sxmgmx = gAlice->Field()->Max();
519   
520   // For polyethilene (CH2) 
521   Float_t ape[2] = { 12., 1. };
522   Float_t zpe[2] = {  6., 1. };
523   Float_t wpe[2] = {  1., 2. };
524   Float_t dpe    = 0.95;
525
526   // For mylar (C5H4O2) 
527   Float_t amy[3] = { 12., 1., 16. };
528   Float_t zmy[3] = {  6., 1.,  8. };
529   Float_t wmy[3] = {  5., 4.,  2. };
530   Float_t dmy    = 1.39;
531
532   // For CO2 
533   Float_t aco[2] = { 12., 16. };
534   Float_t zco[2] = {  6.,  8. };
535   Float_t wco[2] = {  1.,  2. };
536   Float_t dco    = 0.001977;
537
538   // For water
539   Float_t awa[2] = {  1., 16. };
540   Float_t zwa[2] = {  1.,  8. };
541   Float_t wwa[2] = {  2.,  1. };
542   Float_t dwa    = 1.0;
543
544   // For isobutane (C4H10)
545   Float_t ais[2] = { 12.,  1. };
546   Float_t zis[2] = {  6.,  1. };
547   Float_t wis[2] = {  4., 10. };
548   Float_t dis    = 0.00267;
549
550   // For plexiglas (C5H8O2)
551   Float_t apg[3] = { 12.011 ,  1.0    , 15.9994 };
552   Float_t zpg[3] = {  6.0   ,  1.0    ,  8.0    };
553   Float_t wpg[3] = {  5.0   ,  8.0    ,  2.0    };
554   Float_t dpg    = 1.18; 
555
556   // For Xe/CO2-gas-mixture 
557   // Xe-content of the Xe/CO2-mixture (85% / 15%) 
558   Float_t fxc    = .85;
559   // Xe-content of the Xe/Isobutane-mixture (97% / 3%) 
560   Float_t fxi    = .97;
561   Float_t dxe    = .005858;
562   
563   // General tracking parameter
564   Float_t tmaxfd = -10.;
565   Float_t stemax = -1e10;
566   Float_t deemax = -0.1;
567   Float_t epsil  =  1e-4;
568   Float_t stmin  = -0.001;
569   
570   Float_t absl, radl, d, buf[1];
571   Float_t agm[2], zgm[2], wgm[2];
572   Float_t dgm1, dgm2;
573   Int_t   nbuf;
574   
575   //////////////////////////////////////////////////////////////////////////
576   //     Define Materials 
577   //////////////////////////////////////////////////////////////////////////
578
579   AliMaterial( 1, "Al"   ,  26.98, 13.0, 2.7     ,     8.9 ,    37.2);
580   AliMaterial( 2, "Air"  ,  14.61,  7.3, 0.001205, 30420.0 , 67500.0);
581   AliMaterial( 4, "Xe"   , 131.29, 54.0, dxe     ,  1447.59,     0.0);
582   AliMaterial( 5, "Cu"   ,  63.54, 29.0, 8.96    ,     1.43,    14.8);
583   AliMaterial( 6, "C"    ,  12.01,  6.0, 2.265   ,    18.8 ,    74.4);
584   AliMaterial(12, "G10"  ,  20.00, 10.0, 1.7     ,    19.4 ,   999.0);
585   AliMaterial(15, "Sn"   , 118.71, 50.0, 7.31    ,     1.21,    14.8);
586   AliMaterial(16, "Si"   ,  28.09, 14.0, 2.33    ,     9.36,    37.2);
587   AliMaterial(17, "Epoxy",  17.75,  8.9, 1.8     ,    21.82,   999.0);
588
589   // Mixtures 
590   AliMixture(3, "Polyethilene",   ape, zpe, dpe, -2, wpe);
591   AliMixture(7, "Mylar",          amy, zmy, dmy, -3, wmy);
592   AliMixture(8, "CO2",            aco, zco, dco, -2, wco);
593   AliMixture(9, "Isobutane",      ais, zis, dis, -2, wis);
594   AliMixture(13,"Water",          awa, zwa, dwa, -2, wwa);
595   AliMixture(14,"Plexiglas",      apg, zpg, dpg, -3, wpg);
596
597   // Gas mixtures
598   Char_t namate[21]="";
599   // Xe/CO2-mixture
600   // Get properties of Xe 
601   gMC->Gfmate((*fIdmate)[4], namate, agm[0], zgm[0], d, radl, absl, buf, nbuf);
602   // Get properties of CO2 
603   gMC->Gfmate((*fIdmate)[8], namate, agm[1], zgm[1], d, radl, absl, buf, nbuf);
604   // Create gas mixture 
605   wgm[0] = fxc;
606   wgm[1] = 1. - fxc;
607   dgm1   = wgm[0] * dxe + wgm[1] * dco;
608   AliMixture(10, "Gas mixture 1", agm, zgm, dgm1,  2, wgm);
609   // Xe/Isobutane-mixture
610   // Get properties of Xe 
611   gMC->Gfmate((*fIdmate)[4], namate, agm[0], zgm[0], d, radl, absl, buf, nbuf);
612   // Get properties of Isobutane
613   gMC->Gfmate((*fIdmate)[9], namate, agm[1], zgm[1], d, radl, absl, buf, nbuf);
614   // Create gas mixture 
615   wgm[0] = fxi;
616   wgm[1] = 1. - fxi;
617   dgm2   = wgm[0] * dxe + wgm[1] * dis;
618   AliMixture(11, "Gas mixture 2", agm, zgm, dgm2,  2, wgm);
619  
620   //////////////////////////////////////////////////////////////////////////
621   //     Tracking Media Parameters 
622   //////////////////////////////////////////////////////////////////////////
623
624   // Al Frame 
625   AliMedium(1, "Al Frame",   1, 0, isxfld, sxmgmx
626                 , tmaxfd, stemax, deemax, epsil, stmin);
627   // Air 
628   AliMedium(2, "Air",        2, 0, isxfld, sxmgmx
629                 , tmaxfd, stemax, deemax, epsil, stmin);
630   // Polyethilene 
631   AliMedium(3, "Radiator",   3, 0, isxfld, sxmgmx
632                 , tmaxfd, stemax, deemax, epsil, stmin);
633   // Xe 
634   AliMedium(4, "Xe",         4, 1, isxfld, sxmgmx
635                 , tmaxfd, stemax, deemax, epsil, stmin);
636   // Cu pads 
637   AliMedium(5, "Padplane",   5, 1, isxfld, sxmgmx
638                 , tmaxfd, stemax, deemax, epsil, stmin);
639   // Fee + cables 
640   AliMedium(6, "Readout",    1, 0, isxfld, sxmgmx
641                 , tmaxfd, stemax, deemax, epsil, stmin);
642   // C frame 
643   AliMedium(7, "C Frame",    6, 0, isxfld, sxmgmx
644                 , tmaxfd, stemax, deemax, epsil, stmin);
645   // Mylar foils 
646   AliMedium(8, "Mylar",      7, 0, isxfld, sxmgmx
647                 , tmaxfd, stemax, deemax, epsil, stmin);
648   if (fGasMix == 1) {
649     // Gas-mixture (Xe/CO2) 
650     AliMedium(9, "Gas-mix",   10, 1, isxfld, sxmgmx
651                   , tmaxfd, stemax, deemax, epsil, stmin);
652   }
653   else {
654     // Gas-mixture (Xe/Isobutane) 
655     AliMedium(9, "Gas-mix",   11, 1, isxfld, sxmgmx
656                   , tmaxfd, stemax, deemax, epsil, stmin);
657   }
658   // Nomex-honeycomb (use carbon for the time being) 
659   AliMedium(10, "Nomex",      6, 0, isxfld, sxmgmx
660                 , tmaxfd, stemax, deemax, epsil, stmin);
661   // Kapton foils (use Mylar for the time being) 
662   AliMedium(11, "Kapton",     7, 0, isxfld, sxmgmx
663                 , tmaxfd, stemax, deemax, epsil, stmin);
664   // Gas-filling of the radiator 
665   AliMedium(12, "CO2",        8, 0, isxfld, sxmgmx
666                 , tmaxfd, stemax, deemax, epsil, stmin);
667   // G10-plates
668   AliMedium(13, "G10-plates",12, 0, isxfld, sxmgmx
669                 , tmaxfd, stemax, deemax, epsil, stmin);
670   // Cooling water
671   AliMedium(14, "Water",     13, 0, isxfld, sxmgmx
672                 , tmaxfd, stemax, deemax, epsil, stmin);
673   // Rohacell (plexiglas) for the radiator
674   AliMedium(15, "Rohacell",  14, 0, isxfld, sxmgmx
675                 , tmaxfd, stemax, deemax, epsil, stmin);
676   // Al layer in MCMs
677   AliMedium(16, "MCM-Al"  ,   1, 0, isxfld, sxmgmx
678                 , tmaxfd, stemax, deemax, epsil, stmin);
679   // Sn layer in MCMs
680   AliMedium(17, "MCM-Sn"  ,  15, 0, isxfld, sxmgmx
681                 , tmaxfd, stemax, deemax, epsil, stmin);
682   // Cu layer in MCMs
683   AliMedium(18, "MCM-Cu"  ,   5, 0, isxfld, sxmgmx
684                 , tmaxfd, stemax, deemax, epsil, stmin);
685   // G10 layer in MCMs
686   AliMedium(19, "MCM-G10" ,  12, 0, isxfld, sxmgmx
687                 , tmaxfd, stemax, deemax, epsil, stmin);
688   // Si in readout chips
689   AliMedium(20, "Chip-Si" ,  16, 0, isxfld, sxmgmx
690                 , tmaxfd, stemax, deemax, epsil, stmin);
691   // Epoxy in readout chips
692   AliMedium(21, "Chip-Ep" ,  17, 0, isxfld, sxmgmx
693                 , tmaxfd, stemax, deemax, epsil, stmin);
694   // PE in connectors
695   AliMedium(22, "Conn-PE" ,   3, 0, isxfld, sxmgmx
696                 , tmaxfd, stemax, deemax, epsil, stmin);
697   // Cu in connectors
698   AliMedium(23, "Chip-Cu" ,   5, 0, isxfld, sxmgmx
699                 , tmaxfd, stemax, deemax, epsil, stmin);
700   // Al of cooling pipes
701   AliMedium(24, "Cooling" ,   1, 0, isxfld, sxmgmx
702                 , tmaxfd, stemax, deemax, epsil, stmin);
703   // Cu in services
704   AliMedium(25, "Serv-Cu" ,   5, 0, isxfld, sxmgmx
705                 , tmaxfd, stemax, deemax, epsil, stmin);
706
707   // Save the density values for the TRD absorbtion
708   fFoilDensity = dmy;
709   if (fGasMix == 1)
710     fGasDensity = dgm1;
711   else
712     fGasDensity = dgm2;
713
714 }
715
716 //_____________________________________________________________________________
717 void AliTRD::DrawModule() const
718 {
719   //
720   // Draw a shaded view of the Transition Radiation Detector version 0
721   //
722
723   // Set everything unseen
724   gMC->Gsatt("*"   ,"SEEN",-1);
725   
726   // Set ALIC mother transparent
727   gMC->Gsatt("ALIC","SEEN", 0);
728   
729   // Set the volumes visible
730   if (fGeometry->IsVersion() == 0) {
731     gMC->Gsatt("B071","SEEN", 0);
732     gMC->Gsatt("B074","SEEN", 0);
733     gMC->Gsatt("B075","SEEN", 0);
734     gMC->Gsatt("B077","SEEN", 0);
735     gMC->Gsatt("BTR1","SEEN", 0);
736     gMC->Gsatt("BTR2","SEEN", 0);
737     gMC->Gsatt("BTR3","SEEN", 0);
738     gMC->Gsatt("UTR1","SEEN", 0);
739     gMC->Gsatt("UTR2","SEEN", 0);
740     gMC->Gsatt("UTR3","SEEN", 0);
741   }
742   else {
743     gMC->Gsatt("B071","SEEN", 0);
744     gMC->Gsatt("B074","SEEN", 0);
745     gMC->Gsatt("B075","SEEN", 0);
746     gMC->Gsatt("B077","SEEN", 0);
747     gMC->Gsatt("BTR1","SEEN", 0);
748     gMC->Gsatt("BTR2","SEEN", 0);
749     gMC->Gsatt("BTR3","SEEN", 0);
750     gMC->Gsatt("UTR1","SEEN", 0);
751     if (fGeometry->GetPHOShole())
752       gMC->Gsatt("UTR2","SEEN", 0);
753     if (fGeometry->GetRICHhole())
754       gMC->Gsatt("UTR3","SEEN", 0);
755   }
756 //   gMC->Gsatt("UCII","SEEN", 0);
757 //   gMC->Gsatt("UCIM","SEEN", 0);
758 //   gMC->Gsatt("UCIO","SEEN", 0);
759 //   gMC->Gsatt("UL02","SEEN", 1);
760 //   gMC->Gsatt("UL05","SEEN", 1);
761 //   gMC->Gsatt("UL06","SEEN", 1);
762   
763   gMC->Gdopt("hide", "on");
764   gMC->Gdopt("shad", "on");
765   gMC->Gsatt("*", "fill", 7);
766   gMC->SetClipBox(".");
767   gMC->SetClipBox("*", 0, 2000, -2000, 2000, -2000, 2000);
768   gMC->DefaultRange();
769   gMC->Gdraw("alic", 40, 30, 0, 12, 9.4, .021, .021);
770   gMC->Gdhead(1111, "Transition Radiation Detector");
771   gMC->Gdman(18, 4, "MAN");
772
773 }
774
775 //_____________________________________________________________________________
776 Int_t AliTRD::DistancetoPrimitive(Int_t , Int_t ) const
777 {
778   //
779   // Distance between the mouse and the TRD detector on the screen
780   // Dummy routine
781   
782   return 9999;
783
784 }
785  
786 //_____________________________________________________________________________
787 void AliTRD::Init()
788 {
789   //
790   // Initialize the TRD detector after the geometry has been created
791   //
792
793   Int_t i;
794
795   if (fDebug) {
796     printf("\n%s: ",ClassName());
797     for (i = 0; i < 35; i++) printf("*");
798     printf(" TRD_INIT ");
799     for (i = 0; i < 35; i++) printf("*");
800     printf("\n");
801   }
802
803   if      (fGeometry->IsVersion() == 0) {
804     printf("%s: Geometry for spaceframe with holes initialized\n",ClassName());
805   }
806   else if (fGeometry->IsVersion() == 1) {
807     printf("%s: Geometry for spaceframe without holes initialized\n",ClassName());
808     if (fGeometry->GetPHOShole())
809       printf("%s: Leave space in front of PHOS free\n",ClassName());
810     if (fGeometry->GetRICHhole())
811       printf("%s: Leave space in front of RICH free\n",ClassName());
812   }
813   
814   if (fGasMix == 1) {
815     printf("%s: Gas Mixture: 85%% Xe + 15%% CO2\n",ClassName());
816   }
817   else {
818     printf("%s: Gas Mixture: 97%% Xe + 3%% Isobutane\n",ClassName());
819   }
820
821 }
822
823 //_____________________________________________________________________________
824 void AliTRD::LoadPoints(Int_t )
825 {
826   //
827   // Store x, y, z of all hits in memory.
828   // Hit originating from TR photons are given a different color
829   //
830
831   //if (!fDrawTR) {
832   //  AliDetector::LoadPoints(track);
833   //  return;
834   //}
835
836   if ((fHits == 0) && (fTrackHits == 0)) return;
837
838   Int_t nhits;
839   if (fHitType < 2) {
840     nhits = fHits->GetEntriesFast();
841   }
842   else {
843     nhits = fTrackHits->GetEntriesFast();
844   } 
845   if (nhits == 0) return;
846
847   Int_t tracks = gAlice->GetNtrack();
848   if (fPoints == 0) fPoints = new TObjArray(tracks);
849
850   AliTRDhit *ahit;
851   
852   Int_t    *ntrkE = new Int_t[tracks];
853   Int_t    *ntrkT = new Int_t[tracks];
854   Int_t    *limiE = new Int_t[tracks];
855   Int_t    *limiT = new Int_t[tracks];
856   Float_t **coorE = new Float_t*[tracks];
857   Float_t **coorT = new Float_t*[tracks];
858   for(Int_t i = 0; i < tracks; i++) {
859     ntrkE[i] = 0;
860     ntrkT[i] = 0;
861     coorE[i] = 0;
862     coorT[i] = 0;
863     limiE[i] = 0;
864     limiT[i] = 0;
865   }
866   
867   AliTRDpoints *points = 0;
868   Float_t      *fp     = 0;
869   Int_t         trk;
870   Int_t         chunk  = nhits / 4 + 1;
871
872   // Loop over all the hits and store their position
873   ahit = (AliTRDhit *) FirstHit(-1);
874   while (ahit) {
875
876     // dEdx hits
877     if (ahit->GetCharge() >= 0) {
878
879       trk = ahit->GetTrack();
880       if (ntrkE[trk] == limiE[trk]) {
881         // Initialise a new track
882         fp = new Float_t[3*(limiE[trk]+chunk)];
883         if (coorE[trk]) {
884           memcpy(fp,coorE[trk],sizeof(Float_t)*3*limiE[trk]);
885           delete [] coorE[trk];
886         }
887         limiE[trk] += chunk;
888         coorE[trk]  = fp;
889       } 
890       else {
891         fp = coorE[trk];
892       }
893       fp[3*ntrkE[trk]  ] = ahit->X();
894       fp[3*ntrkE[trk]+1] = ahit->Y();
895       fp[3*ntrkE[trk]+2] = ahit->Z();
896       ntrkE[trk]++;
897
898     }
899     // TR photon hits
900     else if ((ahit->GetCharge() < 0) && (fDrawTR)) {
901
902       trk = ahit->GetTrack();
903       if (ntrkT[trk] == limiT[trk]) {
904         // Initialise a new track
905         fp = new Float_t[3*(limiT[trk]+chunk)];
906         if (coorT[trk]) {
907           memcpy(fp,coorT[trk],sizeof(Float_t)*3*limiT[trk]);
908           delete [] coorT[trk];
909         }
910         limiT[trk] += chunk;
911         coorT[trk]  = fp;
912       } 
913       else {
914         fp = coorT[trk];
915       }
916       fp[3*ntrkT[trk]  ] = ahit->X();
917       fp[3*ntrkT[trk]+1] = ahit->Y();
918       fp[3*ntrkT[trk]+2] = ahit->Z();
919       ntrkT[trk]++;
920
921     }
922
923     ahit = (AliTRDhit *) NextHit();
924
925   }
926
927   for (trk = 0; trk < tracks; ++trk) {
928
929     if (ntrkE[trk] || ntrkT[trk]) {
930
931       points = new AliTRDpoints();
932       points->SetDetector(this);
933       points->SetParticle(trk);
934
935       // Set the dEdx points
936       if (ntrkE[trk]) {
937         points->SetMarkerColor(GetMarkerColor());
938         points->SetMarkerSize(GetMarkerSize());
939         points->SetPolyMarker(ntrkE[trk],coorE[trk],GetMarkerStyle());
940         delete [] coorE[trk];
941         coorE[trk] = 0;
942       }
943
944       // Set the TR photon points
945       if (ntrkT[trk]) {
946         points->SetTRpoints(ntrkT[trk],coorT[trk]);
947         delete [] coorT[trk];
948         coorT[trk] = 0;
949       }
950
951       fPoints->AddAt(points,trk);
952
953     }
954
955   }
956
957   delete [] coorE;
958   delete [] coorT;
959   delete [] ntrkE;
960   delete [] ntrkT;
961   delete [] limiE;
962   delete [] limiT;
963
964 }
965
966 //_____________________________________________________________________________
967 void AliTRD::MakeBranch(Option_t* option)
968 {
969   //
970   // Create Tree branches for the TRD digits.
971   //
972
973   Int_t  buffersize = 4000;
974   Char_t branchname[15];
975   sprintf(branchname,"%s",GetName());
976
977   const char *cD = strstr(option,"D");
978
979   AliDetector::MakeBranch(option);
980
981   if (fDigits && gAlice->TreeD() && cD) {
982     MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,buffersize,0);
983   }     
984
985   if (fHitType > 1) {
986     MakeBranch2(option,0); 
987   }
988
989 }
990
991 //_____________________________________________________________________________
992 void AliTRD::ResetDigits()
993 {
994   //
995   // Reset number of digits and the digits array for this detector
996   //
997
998   fNdigits = 0;
999   if (fDigits) fDigits->Clear();
1000
1001 }
1002
1003 //_____________________________________________________________________________
1004 void AliTRD::ResetRecPoints()
1005 {
1006   //
1007   // Reset number of reconstructed points and the point array
1008   //
1009
1010   if (fRecPoints) {
1011     fNRecPoints = 0;
1012     Int_t nentr = fRecPoints->GetEntriesFast();
1013     for (Int_t i = 0; i < nentr; i++) delete fRecPoints->RemoveAt(i);
1014   }
1015
1016 }
1017
1018 //_____________________________________________________________________________
1019 void AliTRD::SetTreeAddress()
1020 {
1021   //
1022   // Set the branch addresses for the trees.
1023   //
1024
1025   Char_t branchname[15];
1026
1027   if ( fLoader->TreeH() && (fHits == 0x0))  fHits = new TClonesArray("AliTRDhit",405);
1028   AliDetector::SetTreeAddress();
1029
1030   TBranch *branch;
1031   TTree   *treeR = fLoader->TreeR();
1032
1033   if (treeR) {
1034     sprintf(branchname,"%scluster",GetName());
1035     if (fRecPoints == 0x0) fRecPoints  = new TObjArray(400);
1036     if (fRecPoints) {
1037       branch = treeR->GetBranch(branchname);
1038       if (branch) {
1039         branch->SetAddress(&fRecPoints);
1040       }
1041     }
1042   }
1043
1044   if (fHitType > 0) {
1045     SetTreeAddress2();    
1046   }
1047
1048 }
1049
1050 //_____________________________________________________________________________
1051 void AliTRD::SetGasMix(Int_t imix)
1052 {
1053   //
1054   // Defines the gas mixture (imix=0:  Xe/Isobutane imix=1: Xe/CO2)
1055   //
1056   
1057   if ((imix < 0) || (imix > 1)) {
1058     printf("Wrong input value: %d\n",imix);
1059     printf("Use standard setting\n");
1060     fGasMix = 1;
1061     return;
1062   }
1063
1064   fGasMix = imix;
1065
1066 }
1067
1068 //_____________________________________________________________________________
1069 void AliTRD::SetPHOShole()
1070 {
1071   //
1072   // Selects a geometry with a hole in front of the PHOS
1073   //
1074
1075   fGeometry->SetPHOShole();
1076
1077 }
1078
1079 //_____________________________________________________________________________
1080 void AliTRD::SetRICHhole()
1081 {
1082   //
1083   // Selects a geometry with a hole in front of the RICH
1084   //
1085
1086   fGeometry->SetRICHhole();
1087
1088 }
1089
1090 //_____________________________________________________________________________
1091 AliTRD &AliTRD::operator=(const AliTRD &trd)
1092 {
1093   //
1094   // Assignment operator
1095   //
1096
1097   if (this != &trd) ((AliTRD &) trd).Copy(*this);
1098   return *this;
1099
1100
1101
1102 //_____________________________________________________________________________
1103 void AliTRD::FinishPrimary()
1104 {
1105   //
1106   // Store the hits in the containers after all primaries are finished
1107   //
1108
1109   if (fTrackHits) { 
1110     fTrackHits->FlushHitStack();
1111   }
1112
1113 }
1114
1115 //_____________________________________________________________________________
1116 void AliTRD::RemapTrackHitIDs(Int_t *map)
1117 {
1118   //
1119   // Remap the track IDs
1120   //
1121
1122   if (!fTrackHits) {
1123     return;
1124   }
1125
1126   if (fTrackHits) {
1127     TClonesArray *arr = fTrackHits->GetArray();;
1128     for (Int_t i = 0; i < arr->GetEntriesFast(); i++){
1129       AliTrackHitsParamV2 *info = (AliTrackHitsParamV2 *) (arr->At(i));
1130       info->fTrackID = map[info->fTrackID];
1131     }
1132   }
1133
1134 }
1135
1136 //_____________________________________________________________________________
1137 void AliTRD::ResetHits()
1138 {
1139   //
1140   // Reset the hits
1141   //
1142
1143   AliDetector::ResetHits();
1144   if (fTrackHits) {
1145     fTrackHits->Clear();
1146   }
1147
1148 }
1149
1150 //_____________________________________________________________________________
1151 AliHit* AliTRD::FirstHit(Int_t track)
1152 {
1153   //
1154   // Return the first hit of a track
1155   //
1156
1157   if (fHitType > 1) {
1158     return FirstHit2(track);
1159   }
1160
1161   return AliDetector::FirstHit(track);
1162
1163 }
1164
1165 //_____________________________________________________________________________
1166 AliHit* AliTRD::NextHit()
1167 {
1168   //
1169   // Returns the next hit of a track
1170   //
1171
1172   if (fHitType > 1) {
1173     return NextHit2();
1174   }
1175
1176   return AliDetector::NextHit();
1177
1178 }
1179
1180 //_____________________________________________________________________________
1181 AliHit* AliTRD::FirstHit2(Int_t track) 
1182 {
1183   //
1184   // Initializes the hit iterator.
1185   // Returns the address of the first hit of a track.
1186   // If <track> >= 0 the track is read from disk,
1187   // while if <track> < 0 the first hit of the current
1188   // track is returned.
1189   //
1190
1191   if (track >= 0) {
1192     gAlice->ResetHits();
1193     TreeH()->GetEvent(track);
1194   }
1195   
1196   if (fTrackHits) {
1197     fTrackHits->First();
1198     return (AliHit*) fTrackHits->GetHit();
1199   }
1200   else {
1201     return 0;
1202   }
1203
1204 }
1205
1206 //_____________________________________________________________________________
1207 AliHit* AliTRD::NextHit2()
1208 {
1209   //
1210   // Returns the next hit of the current track
1211   //
1212
1213   if (fTrackHits) {
1214     fTrackHits->Next();
1215     return (AliHit *) fTrackHits->GetHit();
1216   }
1217   else {
1218     return 0;
1219   }
1220
1221 }
1222
1223 //_____________________________________________________________________________
1224 void AliTRD::MakeBranch2(Option_t *option, const char* )
1225 {
1226   //
1227   // Create a new branch in the current Root tree.
1228   // The branch of fHits is automatically split.
1229   //
1230
1231   if (fHitType < 2) {
1232     return;
1233   }
1234
1235   char branchname[10];
1236   sprintf(branchname,"%s2",GetName());
1237
1238   // Get the pointer to the header
1239   const char *cH = strstr(option,"H");
1240  
1241   if (!fTrackHits) {
1242     fTrackHits = new AliTRDtrackHits();
1243    }
1244
1245
1246   if (fTrackHits && TreeH() && cH) 
1247    {
1248     TreeH()->Branch(branchname,"AliTRDtrackHits",&fTrackHits,fBufferSize,99);
1249     Info("MakeBranch2","Making Branch %s for trackhits",branchname);
1250    }
1251 }
1252
1253 //_____________________________________________________________________________
1254 void AliTRD::SetTreeAddress2()
1255 {
1256   //
1257   // Set the branch address for the trackHits tree
1258   //
1259   if (GetDebug())  Info("SetTreeAddress2","");
1260
1261   TBranch *branch;
1262   char branchname[20];
1263   sprintf(branchname,"%s2",GetName());
1264   
1265   // Branch address for hit tree
1266   TTree *treeH = TreeH();
1267   if ((treeH) && (fHitType > 0)) {
1268     branch = treeH->GetBranch(branchname);
1269     if (branch) 
1270      {
1271        branch->SetAddress(&fTrackHits);
1272        if (GetDebug())
1273          Info("SetTreeAddress2","Success.");
1274      }
1275     else
1276      {
1277        if (GetDebug())
1278          Info("SetTreeAddress2","Can NOT get the branch %s",branchname);
1279      }
1280   }
1281
1282 }
1283
1284 //_____________________________________________________________________________
1285 void AliTRD::AddHit2(Int_t track, Int_t det, Float_t *hits, Int_t q
1286                    , Bool_t inDrift)
1287 {
1288   //
1289   // Add a hit to the list
1290   //
1291
1292   Int_t rtrack;
1293
1294   if (fIshunt) {
1295     Int_t primary = gAlice->GetPrimary(track);
1296     gAlice->Particle(primary)->SetBit(kKeepBit);
1297     rtrack = primary;
1298   } 
1299   else {
1300     rtrack = track;
1301     gAlice->FlagTrack(track);
1302   }
1303
1304   if ((fTrackHits) && (fHitType > 0)) {
1305     fTrackHits->AddHitTRD(det,rtrack,hits[0],hits[1],hits[2],q,inDrift);
1306   }
1307
1308 }
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647