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