]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRD.cxx
Removing AliMC and AliMCProcess
[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 /*
17 $Log$
18 Revision 1.42  2002/10/22 15:53:08  alibrary
19 Introducing Riostream.h
20
21 Revision 1.41  2002/10/14 14:57:43  hristov
22 Merging the VirtualMC branch to the main development branch (HEAD)
23
24 Revision 1.36.6.2  2002/07/24 10:09:30  alibrary
25 Updating VirtualMC
26
27 Revision 1.40  2002/06/13 08:11:56  cblume
28 Add the track references
29
30 Revision 1.39  2002/06/12 09:54:35  cblume
31 Update of tracking code provided by Sergei
32
33 Revision 1.38  2002/03/28 14:59:07  cblume
34 Coding conventions
35
36 Revision 1.37  2002/03/25 20:01:49  cblume
37 Introduce parameter class
38
39 Revision 1.36  2002/02/11 14:25:27  cblume
40 Geometry update, compressed hit structure
41
42 Revision 1.35  2001/11/14 12:08:44  cblume
43 Remove unneccessary header files
44
45 Revision 1.34  2001/11/14 10:50:45  cblume
46 Changes in digits IO. Add merging of summable digits
47
48 Revision 1.33  2001/11/06 17:19:41  cblume
49 Add detailed geometry and simple simulator
50
51 Revision 1.32  2001/10/08 06:57:33  hristov
52 Branches for  TRD digits are created only during the digitisation
53
54 Revision 1.31  2001/08/30 09:30:30  hristov
55 The split level of branches is set to 99
56
57 Revision 1.30  2001/05/28 17:07:58  hristov
58 Last minute changes; ExB correction in AliTRDclusterizerV1; taking into account of material in G10 TEC frames and material between TEC planes (C.Blume,S.Sedykh)
59
60 Revision 1.29  2001/05/21 16:45:47  hristov
61 Last minute changes (C.Blume)
62
63 Revision 1.28  2001/05/16 14:57:27  alibrary
64 New files for folders and Stack
65
66 Revision 1.27  2001/05/08 07:05:02  hristov
67 Loop variable declared once (HP, Sun)
68
69 Revision 1.26  2001/05/07 08:03:22  cblume
70 Generate also hits in the amplification region
71
72 Revision 1.25  2001/03/13 09:30:35  cblume
73 Update of digitization. Moved digit branch definition to AliTRD
74
75 Revision 1.24  2001/01/26 19:56:49  hristov
76 Major upgrade of AliRoot code
77
78 Revision 1.23  2000/11/01 14:53:20  cblume
79 Merge with TRD-develop
80
81 Revision 1.17.2.6  2000/10/15 23:29:08  cblume
82 Introduced more detailed geometry for the display
83
84 Revision 1.17.2.5  2000/10/06 16:49:46  cblume
85 Made Getters const
86
87 Revision 1.17.2.4  2000/10/04 16:34:57  cblume
88 Replace include files by forward declarations
89
90 Revision 1.17.2.3  2000/09/22 14:45:17  cblume
91 Included changes for the tracking
92
93 Revision 1.17.2.2  2000/09/18 13:25:13  cblume
94 Included LoadPoints() method to display the TR photons
95
96 Revision 1.22  2000/10/02 21:28:19  fca
97 Removal of useless dependecies via forward declarations
98
99 Revision 1.21  2000/06/09 11:10:07  cblume
100 Compiler warnings and coding conventions, next round
101
102 Revision 1.20  2000/06/08 18:32:57  cblume
103 Make code compliant to coding conventions
104
105 Revision 1.19  2000/06/07 16:25:37  cblume
106 Try to remove compiler warnings on Sun and HP
107
108 Revision 1.18  2000/05/08 16:17:27  cblume
109 Merge TRD-develop
110
111 Revision 1.21  2000/06/09 11:10:07  cblume
112 Compiler warnings and coding conventions, next round
113
114 Revision 1.20  2000/06/08 18:32:57  cblume
115 Make code compliant to coding conventions
116
117 Revision 1.19  2000/06/07 16:25:37  cblume
118 Try to remove compiler warnings on Sun and HP
119
120 Revision 1.18  2000/05/08 16:17:27  cblume
121 Merge TRD-develop
122
123 Revision 1.17.2.1  2000/05/08 14:28:59  cblume
124 Introduced SetPHOShole() and SetRICHhole(). AliTRDrecPoint container is now a TObjArray
125
126 Revision 1.17  2000/02/28 19:10:26  cblume
127 Include the new TRD classes
128
129 Revision 1.16.2.2  2000/02/28 17:53:24  cblume
130 Introduce TRD geometry classes
131
132 Revision 1.16.2.1  2000/02/28 17:04:19  cblume
133 Include functions and data members for AliTRDrecPoint
134
135 Revision 1.16  2000/01/19 17:17:35  fca
136 Introducing a list of lists of hits -- more hits allowed for detector now
137
138 Revision 1.15  1999/11/02 17:04:25  fca
139 Small syntax change for HP compiler
140
141 Revision 1.14  1999/11/02 16:57:02  fca
142 Avoid non ansi warnings on HP compilers
143
144 Revision 1.13  1999/11/02 16:35:56  fca
145 New version of TRD introduced
146
147 Revision 1.12  1999/11/01 20:41:51  fca
148 Added protections against using the wrong version of FRAME
149
150 Revision 1.11  1999/09/29 09:24:34  fca
151 Introduction of the Copyright and cvs Log
152
153 */
154
155 ///////////////////////////////////////////////////////////////////////////////
156 //                                                                           //
157 //  Transition Radiation Detector                                            //
158 //  This class contains the basic functions for the Transition Radiation     //
159 //  Detector.                                                                //
160 //                                                                           //
161 ///////////////////////////////////////////////////////////////////////////////
162
163 #include <stdlib.h>
164 #include <Riostream.h>
165
166 #include <TMath.h>
167 #include <TNode.h>
168 #include <TGeometry.h>
169 #include <TTree.h>                                                              
170 #include <TPGON.h> 
171 #include <TFile.h>
172 #include <TROOT.h>
173 #include <TParticle.h>
174 #include <TLorentzVector.h>
175
176 #include "AliRun.h"
177 #include "AliConst.h"
178 #include "AliDigit.h"
179 #include "AliMagF.h"
180
181 #include "AliTrackReference.h"
182  
183 #include "AliTRD.h"
184 #include "AliTRDhit.h"
185 #include "AliTRDpoints.h"
186 #include "AliTRDdigit.h"
187 #include "AliTRDdigitizer.h"
188 #include "AliTRDclusterizer.h"
189 #include "AliTRDgeometryHole.h"
190 #include "AliTRDgeometryFull.h"
191 #include "AliTRDrecPoint.h"
192 #include "AliTRDcluster.h"
193 #include "AliTRDdigitsManager.h"
194 #include "AliTRDtrackHits.h"  
195
196 ClassImp(AliTRD)
197  
198 //_____________________________________________________________________________
199 AliTRD::AliTRD()
200 {
201   //
202   // Default constructor
203   //
204
205   fIshunt        = 0;
206   fGasMix        = 0;
207   fHits          = 0;
208   fDigits        = 0;
209
210   fRecPoints     = 0;
211   fNRecPoints    = 0;
212
213   fGeometry      = 0;
214
215   fGasDensity    = 0;
216   fFoilDensity   = 0;
217
218   fDrawTR        = 0;
219   fDisplayType   = 0;
220  
221   fTrackHits     = 0; 
222   fHitType       = 0; 
223
224 }
225  
226 //_____________________________________________________________________________
227 AliTRD::AliTRD(const char *name, const char *title)
228        : AliDetector(name,title)
229 {
230   //
231   // Standard constructor for the TRD
232   //
233
234   // Check that FRAME is there otherwise we have no place where to
235   // put TRD
236   AliModule* frame = gAlice->GetModule("FRAME");
237   if (!frame) {
238     Error("Ctor","TRD needs FRAME to be present\n");
239     exit(1);
240   } 
241
242   // Define the TRD geometry according to the FRAME geometry
243   if      (frame->IsVersion() == 0) {
244     // Geometry with hole
245     fGeometry = new AliTRDgeometryHole();
246   }
247   else if (frame->IsVersion() == 1) {
248     // Geometry without hole
249     fGeometry = new AliTRDgeometryFull();
250   }
251   else {
252     Error("Ctor","Could not find valid FRAME version\n");
253     exit(1);
254   }
255
256   // Allocate the hit array
257   fHits           = new TClonesArray("AliTRDhit"     ,405);
258   gAlice->AddHitList(fHits);
259
260   // Allocate the digits array
261   fDigits         = 0;
262
263   // Allocate the rec point array
264   fRecPoints     = new TObjArray(400);
265   fNRecPoints    = 0;
266    
267   fIshunt        = 0;
268   fGasMix        = 1;
269
270   fGasDensity    = 0;
271   fFoilDensity   = 0;
272
273   fDrawTR        = 0;
274   fDisplayType   = 0;
275
276   fTrackHits     = 0;
277   fHitType       = 2;
278
279   SetMarkerColor(kWhite);   
280
281 }
282
283 //_____________________________________________________________________________
284 AliTRD::AliTRD(const AliTRD &trd)
285 {
286   //
287   // Copy constructor
288   //
289
290   ((AliTRD &) trd).Copy(*this);
291
292 }
293
294 //_____________________________________________________________________________
295 AliTRD::~AliTRD()
296 {
297   //
298   // TRD destructor
299   //
300
301   fIshunt = 0;
302
303   if (fGeometry) {
304     delete fGeometry;
305     fGeometry  = 0;
306   }
307   if (fHits) {
308     delete fHits;
309     fHits      = 0;
310   }
311   if (fRecPoints) {
312     delete fRecPoints;
313     fRecPoints = 0;
314   }
315   if (fTrackHits) {
316     delete fTrackHits;
317     fTrackHits = 0;
318   }
319
320 }
321
322 //_____________________________________________________________________________
323 void AliTRD::AddCluster(Float_t *pos, Int_t det, Float_t amp
324                       , Int_t *tracks, Float_t *sig, Int_t iType)
325 {
326   //
327   // Add a cluster for the TRD
328   //
329
330   AliTRDcluster *c = new AliTRDcluster();
331
332   c->SetDetector(det);
333   c->AddTrackIndex(tracks);
334   c->SetQ(amp);
335   c->SetY(pos[0]);
336   c->SetZ(pos[1]);
337   c->SetSigmaY2(sig[0]);   
338   c->SetSigmaZ2(sig[1]);
339   c->SetLocalTimeBin(((Int_t) pos[2]));
340
341   switch (iType) {
342   case 0:
343     c->Set2pad();
344     break;
345   case 1:
346     c->Set3pad();
347     break;
348   case 2:
349     c->Set4pad();
350     break;
351   case 3:
352     c->Set5pad();
353     break;
354   case 4:
355     c->SetLarge();
356     break;
357   };
358
359   fRecPoints->Add(c);
360
361 }
362
363 //_____________________________________________________________________________
364 void  AliTRD::AddTrackReference(Int_t label, TLorentzVector p, TLorentzVector x)
365 {
366   //
367   // Add a trackrefernce to the list
368   //
369
370   if (!fTrackReferences) {
371     Error("AddTrackReference","Container fTrackRefernce not active\n");
372     return;
373   }
374
375   Int_t nref = fTrackReferences->GetEntriesFast();
376   TClonesArray &lref = *fTrackReferences;
377   AliTrackReference * ref =  new(lref[nref]) AliTrackReference();
378   ref->SetMomentum(p[0],p[1],p[2]);
379   ref->SetPosition(x[0],x[1],x[2]);
380   ref->SetTrack(label);
381
382 }
383
384 //_____________________________________________________________________________
385 void AliTRD::Hits2Digits()
386 {
387   //
388   // Create digits
389   //
390
391   AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
392                                                   ,"TRD digitizer class");
393   digitizer->SetDebug(GetDebug());
394   digitizer->SetEvent(gAlice->GetEvNumber());
395
396   // Initialization
397   digitizer->InitDetector();
398     
399   // Create the digits
400   digitizer->MakeDigits();
401   
402   // Write the digits into the input file
403   if (digitizer->MakeBranch(fDigitsFile)) {
404
405     digitizer->WriteDigits();
406
407     // Save the digitizer class in the AliROOT 
408     digitizer->Write();
409
410   }
411
412 }
413
414 //_____________________________________________________________________________
415 void AliTRD::Hits2SDigits()
416 {
417   //
418   // Create summable digits
419   //
420
421   AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
422                                                   ,"TRD digitizer class");
423   digitizer->SetDebug(GetDebug());
424
425   // For the summable digits
426   digitizer->SetSDigits(kTRUE);
427   digitizer->SetEvent(gAlice->GetEvNumber());
428
429   // Initialization
430   digitizer->InitDetector();
431     
432   // Create the TRD s-digits branch
433   digitizer->MakeDigits();
434   
435   // Write the digits into the input file
436   if (digitizer->MakeBranch(fDigitsFile)) {
437
438     digitizer->WriteDigits();
439
440     // Save the digitizer class in the AliROOT 
441     digitizer->Write();
442
443   }
444
445 }
446
447 //_____________________________________________________________________________
448 void AliTRD::SDigits2Digits()
449 {
450   //
451   // Create final digits from summable digits
452   //
453
454    // Create the TRD digitizer
455   AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
456                                                   ,"TRD digitizer class");  
457   digitizer->SetDebug(GetDebug());
458
459   // Set the parameter
460   digitizer->SetEvent(gAlice->GetEvNumber());
461
462   // Initialization
463   digitizer->InitDetector();
464
465   // Read the s-digits via digits manager
466   AliTRDdigitsManager *sdigitsManager = new AliTRDdigitsManager();
467   sdigitsManager->SetDebug(GetDebug());
468   sdigitsManager->SetSDigits(kTRUE);
469   if (fDigitsFile) {
470     sdigitsManager->Open(fDigitsFile);
471   }
472   sdigitsManager->CreateArrays();
473   sdigitsManager->ReadDigits();
474
475   // Add the s-digits to the input list 
476   digitizer->AddSDigitsManager(sdigitsManager);
477
478   // Convert the s-digits to normal digits
479   digitizer->SDigits2Digits();
480
481   // Store the digits
482   if (digitizer->MakeBranch(fDigitsFile)) {
483
484     digitizer->WriteDigits();
485
486   }
487
488 }
489
490 //_____________________________________________________________________________
491 void AliTRD::AddHit(Int_t track, Int_t det, Float_t *hits, Int_t q
492                   , Bool_t inDrift)
493 {
494   //
495   // Add a hit for the TRD
496   // 
497   // The data structure is set according to fHitType:
498   //   bit0: standard TClonesArray
499   //   bit1: compressed trackHits structure
500   //
501
502   if (fHitType & 1) {
503     TClonesArray &lhits = *fHits;
504     new(lhits[fNhits++]) AliTRDhit(fIshunt,track,det,hits,q);
505   }
506
507   if (fHitType > 1) {
508     AddHit2(track,det,hits,q,inDrift);
509   }
510
511 }
512
513 //_____________________________________________________________________________
514 void AliTRD::BuildGeometry()
515 {
516   //
517   // Create the ROOT TNode geometry for the TRD
518   //
519
520   TNode *node, *top;
521   TPGON *pgon;
522
523   Float_t rmin, rmax;
524   Float_t zmax1, zmax2;
525
526   Int_t   iPlan;
527  
528   const Int_t kColorTRD = 46;
529   
530   // Find the top node alice
531   top = gAlice->GetGeometry()->GetNode("alice");
532   
533   if      (fDisplayType == 0) {
534
535     pgon = new TPGON("S_TRD","TRD","void",0,360,AliTRDgeometry::Nsect(),4);
536     rmin = AliTRDgeometry::Rmin();
537     rmax = AliTRDgeometry::Rmax();
538     pgon->DefineSection(0,-AliTRDgeometry::Zmax1(),rmax,rmax);
539     pgon->DefineSection(1,-AliTRDgeometry::Zmax2(),rmin,rmax);
540     pgon->DefineSection(2, AliTRDgeometry::Zmax2(),rmin,rmax);
541     pgon->DefineSection(3, AliTRDgeometry::Zmax1(),rmax,rmax);
542     top->cd();
543     node = new TNode("TRD","TRD","S_TRD",0,0,0,"");
544     node->SetLineColor(kColorTRD);
545     fNodes->Add(node);
546
547   }
548   else if (fDisplayType == 1) {
549
550     Char_t name[7];
551
552     Float_t slope = (AliTRDgeometry::Zmax1() - AliTRDgeometry::Zmax2())
553                   / (AliTRDgeometry::Rmax()  - AliTRDgeometry::Rmin());
554
555     rmin  = AliTRDgeometry::Rmin() + AliTRDgeometry::CraHght();
556     rmax  = rmin                   + AliTRDgeometry::CdrHght();
557
558     Float_t thickness = rmin - AliTRDgeometry::Rmin();
559     zmax2 = AliTRDgeometry::Zmax2() + slope * thickness;
560     zmax1 = zmax2 + slope * AliTRDgeometry::DrThick();
561
562     for (iPlan = 0; iPlan < AliTRDgeometry::Nplan(); iPlan++) {
563
564       sprintf(name,"S_TR1%d",iPlan);
565       pgon  = new TPGON(name,"TRD","void",0,360,AliTRDgeometry::Nsect(),4);
566       pgon->DefineSection(0,-zmax1,rmax,rmax);
567       pgon->DefineSection(1,-zmax2,rmin,rmax);
568       pgon->DefineSection(2, zmax2,rmin,rmax);
569       pgon->DefineSection(3, zmax1,rmax,rmax);
570       top->cd();
571       node = new TNode("TRD","TRD",name,0,0,0,"");
572       node->SetLineColor(kColorTRD);
573       fNodes->Add(node);
574
575       Float_t height = AliTRDgeometry::Cheight() + AliTRDgeometry::Cspace(); 
576       rmin  = rmin  + height;
577       rmax  = rmax  + height;
578       zmax1 = zmax1 + slope * height;
579       zmax2 = zmax2 + slope * height;
580
581     }
582
583     thickness += AliTRDgeometry::DrThick();
584     rmin  = AliTRDgeometry::Rmin() + thickness;
585     rmax  = rmin + AliTRDgeometry::AmThick();
586     zmax2 = AliTRDgeometry::Zmax2() + slope * thickness;
587     zmax1 = zmax2 + slope * AliTRDgeometry::AmThick();
588
589     for (iPlan = 0; iPlan < AliTRDgeometry::Nplan(); iPlan++) {
590
591       sprintf(name,"S_TR2%d",iPlan);
592       pgon  = new TPGON(name,"TRD","void",0,360,AliTRDgeometry::Nsect(),4);
593       pgon->DefineSection(0,-zmax1,rmax,rmax);
594       pgon->DefineSection(1,-zmax2,rmin,rmax);
595       pgon->DefineSection(2, zmax2,rmin,rmax);
596       pgon->DefineSection(3, zmax1,rmax,rmax);
597       top->cd();
598       node = new TNode("TRD","TRD",name,0,0,0,"");
599       node->SetLineColor(kColorTRD);
600       fNodes->Add(node);
601
602       Float_t height = AliTRDgeometry::Cheight() + AliTRDgeometry::Cspace(); 
603       rmin  = rmin  + height;
604       rmax  = rmax  + height;
605       zmax1 = zmax1 + slope * height;
606       zmax2 = zmax2 + slope * height;
607
608     }
609
610   }
611
612 }
613  
614 //_____________________________________________________________________________
615 void AliTRD::Copy(TObject &trd)
616 {
617   //
618   // Copy function
619   //
620
621   ((AliTRD &) trd).fGasMix      = fGasMix;
622   ((AliTRD &) trd).fGeometry    = fGeometry;       
623   ((AliTRD &) trd).fRecPoints   = fRecPoints;
624   ((AliTRD &) trd).fNRecPoints  = fNRecPoints;
625   ((AliTRD &) trd).fGasDensity  = fGasDensity;
626   ((AliTRD &) trd).fFoilDensity = fFoilDensity;
627   ((AliTRD &) trd).fDrawTR      = fDrawTR;
628   ((AliTRD &) trd).fDisplayType = fDisplayType;
629   ((AliTRD &) trd).fHitType     = fHitType;
630
631   //AliDetector::Copy(trd);
632
633 }
634
635 //_____________________________________________________________________________
636 void AliTRD::CreateGeometry()
637 {
638   //
639   // Creates the volumes for the TRD chambers
640   //
641
642   // Check that FRAME is there otherwise we have no place where to put the TRD
643   AliModule* frame = gAlice->GetModule("FRAME");
644   if (!frame) {
645     printf(" The TRD needs the FRAME to be defined first\n");
646     return;
647   }
648
649   fGeometry->CreateGeometry(fIdtmed->GetArray() - 1299);
650
651 }
652  
653 //_____________________________________________________________________________
654 void AliTRD::CreateMaterials()
655 {
656   //
657   // Create the materials for the TRD
658   // Origin Y.Foka
659   //
660
661   Int_t   isxfld = gAlice->Field()->Integ();
662   Float_t sxmgmx = gAlice->Field()->Max();
663   
664   // For polyethilene (CH2) 
665   Float_t ape[2] = { 12., 1. };
666   Float_t zpe[2] = {  6., 1. };
667   Float_t wpe[2] = {  1., 2. };
668   Float_t dpe    = 0.95;
669
670   // For mylar (C5H4O2) 
671   Float_t amy[3] = { 12., 1., 16. };
672   Float_t zmy[3] = {  6., 1.,  8. };
673   Float_t wmy[3] = {  5., 4.,  2. };
674   Float_t dmy    = 1.39;
675
676   // For CO2 
677   Float_t aco[2] = { 12., 16. };
678   Float_t zco[2] = {  6.,  8. };
679   Float_t wco[2] = {  1.,  2. };
680   Float_t dco    = 0.001977;
681
682   // For water
683   Float_t awa[2] = {  1., 16. };
684   Float_t zwa[2] = {  1.,  8. };
685   Float_t wwa[2] = {  2.,  1. };
686   Float_t dwa    = 1.0;
687
688   // For isobutane (C4H10)
689   Float_t ais[2] = { 12.,  1. };
690   Float_t zis[2] = {  6.,  1. };
691   Float_t wis[2] = {  4., 10. };
692   Float_t dis    = 0.00267;
693
694   // For plexiglas (C5H8O2)
695   Float_t apg[3] = { 12.011 ,  1.0    , 15.9994 };
696   Float_t zpg[3] = {  6.0   ,  1.0    ,  8.0    };
697   Float_t wpg[3] = {  5.0   ,  8.0    ,  2.0    };
698   Float_t dpg    = 1.18; 
699
700   // For Xe/CO2-gas-mixture 
701   // Xe-content of the Xe/CO2-mixture (85% / 15%) 
702   Float_t fxc    = .85;
703   // Xe-content of the Xe/Isobutane-mixture (97% / 3%) 
704   Float_t fxi    = .97;
705   Float_t dxe    = .005858;
706   
707   // General tracking parameter
708   Float_t tmaxfd = -10.;
709   Float_t stemax = -1e10;
710   Float_t deemax = -0.1;
711   Float_t epsil  =  1e-4;
712   Float_t stmin  = -0.001;
713   
714   Float_t absl, radl, d, buf[1];
715   Float_t agm[2], zgm[2], wgm[2];
716   Float_t dgm1, dgm2;
717   Int_t   nbuf;
718   
719   //////////////////////////////////////////////////////////////////////////
720   //     Define Materials 
721   //////////////////////////////////////////////////////////////////////////
722
723   AliMaterial( 1, "Al"   ,  26.98, 13.0, 2.7     ,     8.9 ,    37.2);
724   AliMaterial( 2, "Air"  ,  14.61,  7.3, 0.001205, 30420.0 , 67500.0);
725   AliMaterial( 4, "Xe"   , 131.29, 54.0, dxe     ,  1447.59,     0.0);
726   AliMaterial( 5, "Cu"   ,  63.54, 29.0, 8.96    ,     1.43,    14.8);
727   AliMaterial( 6, "C"    ,  12.01,  6.0, 2.265   ,    18.8 ,    74.4);
728   AliMaterial(12, "G10"  ,  20.00, 10.0, 1.7     ,    19.4 ,   999.0);
729   AliMaterial(15, "Sn"   , 118.71, 50.0, 7.31    ,     1.21,    14.8);
730   AliMaterial(16, "Si"   ,  28.09, 14.0, 2.33    ,     9.36,    37.2);
731   AliMaterial(17, "Epoxy",  17.75,  8.9, 1.8     ,    21.82,   999.0);
732
733   // Mixtures 
734   AliMixture(3, "Polyethilene",   ape, zpe, dpe, -2, wpe);
735   AliMixture(7, "Mylar",          amy, zmy, dmy, -3, wmy);
736   AliMixture(8, "CO2",            aco, zco, dco, -2, wco);
737   AliMixture(9, "Isobutane",      ais, zis, dis, -2, wis);
738   AliMixture(13,"Water",          awa, zwa, dwa, -2, wwa);
739   AliMixture(14,"Plexiglas",      apg, zpg, dpg, -3, wpg);
740
741   // Gas mixtures
742   Char_t namate[21];
743   // Xe/CO2-mixture
744   // Get properties of Xe 
745   gMC->Gfmate((*fIdmate)[4], namate, agm[0], zgm[0], d, radl, absl, buf, nbuf);
746   // Get properties of CO2 
747   gMC->Gfmate((*fIdmate)[8], namate, agm[1], zgm[1], d, radl, absl, buf, nbuf);
748   // Create gas mixture 
749   wgm[0] = fxc;
750   wgm[1] = 1. - fxc;
751   dgm1   = wgm[0] * dxe + wgm[1] * dco;
752   AliMixture(10, "Gas mixture 1", agm, zgm, dgm1,  2, wgm);
753   // Xe/Isobutane-mixture
754   // Get properties of Xe 
755   gMC->Gfmate((*fIdmate)[4], namate, agm[0], zgm[0], d, radl, absl, buf, nbuf);
756   // Get properties of Isobutane
757   gMC->Gfmate((*fIdmate)[9], namate, agm[1], zgm[1], d, radl, absl, buf, nbuf);
758   // Create gas mixture 
759   wgm[0] = fxi;
760   wgm[1] = 1. - fxi;
761   dgm2   = wgm[0] * dxe + wgm[1] * dis;
762   AliMixture(11, "Gas mixture 2", agm, zgm, dgm2,  2, wgm);
763  
764   //////////////////////////////////////////////////////////////////////////
765   //     Tracking Media Parameters 
766   //////////////////////////////////////////////////////////////////////////
767
768   // Al Frame 
769   AliMedium(1, "Al Frame",   1, 0, isxfld, sxmgmx
770                 , tmaxfd, stemax, deemax, epsil, stmin);
771   // Air 
772   AliMedium(2, "Air",        2, 0, isxfld, sxmgmx
773                 , tmaxfd, stemax, deemax, epsil, stmin);
774   // Polyethilene 
775   AliMedium(3, "Radiator",   3, 0, isxfld, sxmgmx
776                 , tmaxfd, stemax, deemax, epsil, stmin);
777   // Xe 
778   AliMedium(4, "Xe",         4, 1, isxfld, sxmgmx
779                 , tmaxfd, stemax, deemax, epsil, stmin);
780   // Cu pads 
781   AliMedium(5, "Padplane",   5, 1, isxfld, sxmgmx
782                 , tmaxfd, stemax, deemax, epsil, stmin);
783   // Fee + cables 
784   AliMedium(6, "Readout",    1, 0, isxfld, sxmgmx
785                 , tmaxfd, stemax, deemax, epsil, stmin);
786   // C frame 
787   AliMedium(7, "C Frame",    6, 0, isxfld, sxmgmx
788                 , tmaxfd, stemax, deemax, epsil, stmin);
789   // Mylar foils 
790   AliMedium(8, "Mylar",      7, 0, isxfld, sxmgmx
791                 , tmaxfd, stemax, deemax, epsil, stmin);
792   if (fGasMix == 1) {
793     // Gas-mixture (Xe/CO2) 
794     AliMedium(9, "Gas-mix",   10, 1, isxfld, sxmgmx
795                   , tmaxfd, stemax, deemax, epsil, stmin);
796   }
797   else {
798     // Gas-mixture (Xe/Isobutane) 
799     AliMedium(9, "Gas-mix",   11, 1, isxfld, sxmgmx
800                   , tmaxfd, stemax, deemax, epsil, stmin);
801   }
802   // Nomex-honeycomb (use carbon for the time being) 
803   AliMedium(10, "Nomex",      6, 0, isxfld, sxmgmx
804                 , tmaxfd, stemax, deemax, epsil, stmin);
805   // Kapton foils (use Mylar for the time being) 
806   AliMedium(11, "Kapton",     7, 0, isxfld, sxmgmx
807                 , tmaxfd, stemax, deemax, epsil, stmin);
808   // Gas-filling of the radiator 
809   AliMedium(12, "CO2",        8, 0, isxfld, sxmgmx
810                 , tmaxfd, stemax, deemax, epsil, stmin);
811   // G10-plates
812   AliMedium(13, "G10-plates",12, 0, isxfld, sxmgmx
813                 , tmaxfd, stemax, deemax, epsil, stmin);
814   // Cooling water
815   AliMedium(14, "Water",     13, 0, isxfld, sxmgmx
816                 , tmaxfd, stemax, deemax, epsil, stmin);
817   // Rohacell (plexiglas) for the radiator
818   AliMedium(15, "Rohacell",  14, 0, isxfld, sxmgmx
819                 , tmaxfd, stemax, deemax, epsil, stmin);
820   // Al layer in MCMs
821   AliMedium(16, "MCM-Al"  ,   1, 0, isxfld, sxmgmx
822                 , tmaxfd, stemax, deemax, epsil, stmin);
823   // Sn layer in MCMs
824   AliMedium(17, "MCM-Sn"  ,  15, 0, isxfld, sxmgmx
825                 , tmaxfd, stemax, deemax, epsil, stmin);
826   // Cu layer in MCMs
827   AliMedium(18, "MCM-Cu"  ,   5, 0, isxfld, sxmgmx
828                 , tmaxfd, stemax, deemax, epsil, stmin);
829   // G10 layer in MCMs
830   AliMedium(19, "MCM-G10" ,  12, 0, isxfld, sxmgmx
831                 , tmaxfd, stemax, deemax, epsil, stmin);
832   // Si in readout chips
833   AliMedium(20, "Chip-Si" ,  16, 0, isxfld, sxmgmx
834                 , tmaxfd, stemax, deemax, epsil, stmin);
835   // Epoxy in readout chips
836   AliMedium(21, "Chip-Ep" ,  17, 0, isxfld, sxmgmx
837                 , tmaxfd, stemax, deemax, epsil, stmin);
838   // PE in connectors
839   AliMedium(22, "Conn-PE" ,   3, 0, isxfld, sxmgmx
840                 , tmaxfd, stemax, deemax, epsil, stmin);
841   // Cu in connectors
842   AliMedium(23, "Chip-Cu" ,   5, 0, isxfld, sxmgmx
843                 , tmaxfd, stemax, deemax, epsil, stmin);
844   // Al of cooling pipes
845   AliMedium(24, "Cooling" ,   1, 0, isxfld, sxmgmx
846                 , tmaxfd, stemax, deemax, epsil, stmin);
847
848   // Save the density values for the TRD absorbtion
849   fFoilDensity = dmy;
850   if (fGasMix == 1)
851     fGasDensity = dgm1;
852   else
853     fGasDensity = dgm2;
854
855 }
856
857 //_____________________________________________________________________________
858 void AliTRD::DrawModule() const
859 {
860   //
861   // Draw a shaded view of the Transition Radiation Detector version 0
862   //
863
864   // Set everything unseen
865   gMC->Gsatt("*"   ,"SEEN",-1);
866   
867   // Set ALIC mother transparent
868   gMC->Gsatt("ALIC","SEEN", 0);
869   
870   // Set the volumes visible
871   if (fGeometry->IsVersion() == 0) {
872     gMC->Gsatt("B071","SEEN", 0);
873     gMC->Gsatt("B074","SEEN", 0);
874     gMC->Gsatt("B075","SEEN", 0);
875     gMC->Gsatt("B077","SEEN", 0);
876     gMC->Gsatt("BTR1","SEEN", 0);
877     gMC->Gsatt("BTR2","SEEN", 0);
878     gMC->Gsatt("BTR3","SEEN", 0);
879     gMC->Gsatt("UTR1","SEEN", 0);
880     gMC->Gsatt("UTR2","SEEN", 0);
881     gMC->Gsatt("UTR3","SEEN", 0);
882   }
883   else {
884     gMC->Gsatt("B071","SEEN", 0);
885     gMC->Gsatt("B074","SEEN", 0);
886     gMC->Gsatt("B075","SEEN", 0);
887     gMC->Gsatt("B077","SEEN", 0);
888     gMC->Gsatt("BTR1","SEEN", 0);
889     gMC->Gsatt("BTR2","SEEN", 0);
890     gMC->Gsatt("BTR3","SEEN", 0);
891     gMC->Gsatt("UTR1","SEEN", 0);
892     if (fGeometry->GetPHOShole())
893       gMC->Gsatt("UTR2","SEEN", 0);
894     if (fGeometry->GetRICHhole())
895       gMC->Gsatt("UTR3","SEEN", 0);
896   }
897 //   gMC->Gsatt("UCII","SEEN", 0);
898 //   gMC->Gsatt("UCIM","SEEN", 0);
899 //   gMC->Gsatt("UCIO","SEEN", 0);
900 //   gMC->Gsatt("UL02","SEEN", 1);
901 //   gMC->Gsatt("UL05","SEEN", 1);
902 //   gMC->Gsatt("UL06","SEEN", 1);
903   
904   gMC->Gdopt("hide", "on");
905   gMC->Gdopt("shad", "on");
906   gMC->Gsatt("*", "fill", 7);
907   gMC->SetClipBox(".");
908   gMC->SetClipBox("*", 0, 2000, -2000, 2000, -2000, 2000);
909   gMC->DefaultRange();
910   gMC->Gdraw("alic", 40, 30, 0, 12, 9.4, .021, .021);
911   gMC->Gdhead(1111, "Transition Radiation Detector");
912   gMC->Gdman(18, 4, "MAN");
913
914 }
915
916 //_____________________________________________________________________________
917 Int_t AliTRD::DistancetoPrimitive(Int_t , Int_t ) const
918 {
919   //
920   // Distance between the mouse and the TRD detector on the screen
921   // Dummy routine
922   
923   return 9999;
924
925 }
926  
927 //_____________________________________________________________________________
928 void AliTRD::Init()
929 {
930   //
931   // Initialize the TRD detector after the geometry has been created
932   //
933
934   Int_t i;
935
936   if (fDebug) {
937     printf("\n%s: ",ClassName());
938     for (i = 0; i < 35; i++) printf("*");
939     printf(" TRD_INIT ");
940     for (i = 0; i < 35; i++) printf("*");
941     printf("\n");
942   }
943
944   if      (fGeometry->IsVersion() == 0) {
945     printf("%s: Geometry for spaceframe with holes initialized\n",ClassName());
946   }
947   else if (fGeometry->IsVersion() == 1) {
948     printf("%s: Geometry for spaceframe without holes initialized\n",ClassName());
949     if (fGeometry->GetPHOShole())
950       printf("%s: Leave space in front of PHOS free\n",ClassName());
951     if (fGeometry->GetRICHhole())
952       printf("%s: Leave space in front of RICH free\n",ClassName());
953   }
954   
955   if (fGasMix == 1) {
956     printf("%s: Gas Mixture: 85%% Xe + 15%% CO2\n",ClassName());
957   }
958   else {
959     printf("%s: Gas Mixture: 97%% Xe + 3%% Isobutane\n",ClassName());
960   }
961
962 }
963
964 //_____________________________________________________________________________
965 void AliTRD::LoadPoints(Int_t track)
966 {
967   //
968   // Store x, y, z of all hits in memory.
969   // Hit originating from TR photons are given a different color
970   //
971
972   //if (!fDrawTR) {
973   //  AliDetector::LoadPoints(track);
974   //  return;
975   //}
976
977   if ((fHits == 0) && (fTrackHits == 0)) return;
978
979   Int_t nhits;
980   if (fHitType < 2) {
981     nhits = fHits->GetEntriesFast();
982   }
983   else {
984     nhits = fTrackHits->GetEntriesFast();
985   } 
986   if (nhits == 0) return;
987
988   Int_t tracks = gAlice->GetNtrack();
989   if (fPoints == 0) fPoints = new TObjArray(tracks);
990
991   AliTRDhit *ahit;
992   
993   Int_t    *ntrkE = new Int_t[tracks];
994   Int_t    *ntrkT = new Int_t[tracks];
995   Int_t    *limiE = new Int_t[tracks];
996   Int_t    *limiT = new Int_t[tracks];
997   Float_t **coorE = new Float_t*[tracks];
998   Float_t **coorT = new Float_t*[tracks];
999   for(Int_t i = 0; i < tracks; i++) {
1000     ntrkE[i] = 0;
1001     ntrkT[i] = 0;
1002     coorE[i] = 0;
1003     coorT[i] = 0;
1004     limiE[i] = 0;
1005     limiT[i] = 0;
1006   }
1007   
1008   AliTRDpoints *points = 0;
1009   Float_t      *fp     = 0;
1010   Int_t         trk;
1011   Int_t         chunk  = nhits / 4 + 1;
1012
1013   // Loop over all the hits and store their position
1014   ahit = (AliTRDhit *) FirstHit(-1);
1015   while (ahit) {
1016
1017     // dEdx hits
1018     if (ahit->GetCharge() >= 0) {
1019
1020       trk = ahit->GetTrack();
1021       if (ntrkE[trk] == limiE[trk]) {
1022         // Initialise a new track
1023         fp = new Float_t[3*(limiE[trk]+chunk)];
1024         if (coorE[trk]) {
1025           memcpy(fp,coorE[trk],sizeof(Float_t)*3*limiE[trk]);
1026           delete [] coorE[trk];
1027         }
1028         limiE[trk] += chunk;
1029         coorE[trk]  = fp;
1030       } 
1031       else {
1032         fp = coorE[trk];
1033       }
1034       fp[3*ntrkE[trk]  ] = ahit->X();
1035       fp[3*ntrkE[trk]+1] = ahit->Y();
1036       fp[3*ntrkE[trk]+2] = ahit->Z();
1037       ntrkE[trk]++;
1038
1039     }
1040     // TR photon hits
1041     else if ((ahit->GetCharge() < 0) && (fDrawTR)) {
1042
1043       trk = ahit->GetTrack();
1044       if (ntrkT[trk] == limiT[trk]) {
1045         // Initialise a new track
1046         fp = new Float_t[3*(limiT[trk]+chunk)];
1047         if (coorT[trk]) {
1048           memcpy(fp,coorT[trk],sizeof(Float_t)*3*limiT[trk]);
1049           delete [] coorT[trk];
1050         }
1051         limiT[trk] += chunk;
1052         coorT[trk]  = fp;
1053       } 
1054       else {
1055         fp = coorT[trk];
1056       }
1057       fp[3*ntrkT[trk]  ] = ahit->X();
1058       fp[3*ntrkT[trk]+1] = ahit->Y();
1059       fp[3*ntrkT[trk]+2] = ahit->Z();
1060       ntrkT[trk]++;
1061
1062     }
1063
1064     ahit = (AliTRDhit *) NextHit();
1065
1066   }
1067
1068   for (trk = 0; trk < tracks; ++trk) {
1069
1070     if (ntrkE[trk] || ntrkT[trk]) {
1071
1072       points = new AliTRDpoints();
1073       points->SetDetector(this);
1074       points->SetParticle(trk);
1075
1076       // Set the dEdx points
1077       if (ntrkE[trk]) {
1078         points->SetMarkerColor(GetMarkerColor());
1079         points->SetMarkerSize(GetMarkerSize());
1080         points->SetPolyMarker(ntrkE[trk],coorE[trk],GetMarkerStyle());
1081         delete [] coorE[trk];
1082         coorE[trk] = 0;
1083       }
1084
1085       // Set the TR photon points
1086       if (ntrkT[trk]) {
1087         points->SetTRpoints(ntrkT[trk],coorT[trk]);
1088         delete [] coorT[trk];
1089         coorT[trk] = 0;
1090       }
1091
1092       fPoints->AddAt(points,trk);
1093
1094     }
1095
1096   }
1097
1098   delete [] coorE;
1099   delete [] coorT;
1100   delete [] ntrkE;
1101   delete [] ntrkT;
1102   delete [] limiE;
1103   delete [] limiT;
1104
1105 }
1106
1107 //_____________________________________________________________________________
1108 void AliTRD::MakeBranch(Option_t* option, const char *file)
1109 {
1110   //
1111   // Create Tree branches for the TRD digits.
1112   //
1113
1114   Int_t  buffersize = 4000;
1115   Char_t branchname[15];
1116   sprintf(branchname,"%s",GetName());
1117
1118   const char *cD = strstr(option,"D");
1119
1120   AliDetector::MakeBranch(option,file);
1121
1122   if (fDigits && gAlice->TreeD() && cD) {
1123     MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,buffersize,file);
1124   }     
1125
1126   if (fHitType > 1) {
1127     MakeBranch2(option,file); 
1128   }
1129
1130 }
1131
1132 //_____________________________________________________________________________
1133 void AliTRD::ResetDigits()
1134 {
1135   //
1136   // Reset number of digits and the digits array for this detector
1137   //
1138
1139   fNdigits = 0;
1140   if (fDigits) fDigits->Clear();
1141
1142 }
1143
1144 //_____________________________________________________________________________
1145 void AliTRD::ResetRecPoints()
1146 {
1147   //
1148   // Reset number of reconstructed points and the point array
1149   //
1150
1151   if (fRecPoints) {
1152     fNRecPoints = 0;
1153     Int_t nentr = fRecPoints->GetEntriesFast();
1154     for (Int_t i = 0; i < nentr; i++) delete fRecPoints->RemoveAt(i);
1155   }
1156
1157 }
1158
1159 //_____________________________________________________________________________
1160 void AliTRD::SetTreeAddress()
1161 {
1162   //
1163   // Set the branch addresses for the trees.
1164   //
1165
1166   Char_t branchname[15];
1167
1168   AliDetector::SetTreeAddress();
1169
1170   TBranch *branch;
1171   TTree   *treeR = gAlice->TreeR();
1172
1173   if (treeR) {
1174     sprintf(branchname,"%scluster",GetName());
1175     if (fRecPoints) {
1176       branch = treeR->GetBranch(branchname);
1177       if (branch) {
1178         branch->SetAddress(&fRecPoints);
1179       }
1180     }
1181   }
1182
1183   if (fHitType > 0) {
1184     SetTreeAddress2();    
1185   }
1186
1187 }
1188
1189 //_____________________________________________________________________________
1190 void AliTRD::SetGasMix(Int_t imix)
1191 {
1192   //
1193   // Defines the gas mixture (imix=0:  Xe/Isobutane imix=1: Xe/CO2)
1194   //
1195   
1196   if ((imix < 0) || (imix > 1)) {
1197     printf("Wrong input value: %d\n",imix);
1198     printf("Use standard setting\n");
1199     fGasMix = 1;
1200     return;
1201   }
1202
1203   fGasMix = imix;
1204
1205 }
1206
1207 //_____________________________________________________________________________
1208 void AliTRD::SetPHOShole()
1209 {
1210   //
1211   // Selects a geometry with a hole in front of the PHOS
1212   //
1213
1214   fGeometry->SetPHOShole();
1215
1216 }
1217
1218 //_____________________________________________________________________________
1219 void AliTRD::SetRICHhole()
1220 {
1221   //
1222   // Selects a geometry with a hole in front of the RICH
1223   //
1224
1225   fGeometry->SetRICHhole();
1226
1227 }
1228
1229 //_____________________________________________________________________________
1230 AliTRD &AliTRD::operator=(const AliTRD &trd)
1231 {
1232   //
1233   // Assignment operator
1234   //
1235
1236   if (this != &trd) ((AliTRD &) trd).Copy(*this);
1237   return *this;
1238
1239
1240
1241 //_____________________________________________________________________________
1242 void AliTRD::FinishPrimary()
1243 {
1244   //
1245   // Store the hits in the containers after all primaries are finished
1246   //
1247
1248   if (fTrackHits) { 
1249     fTrackHits->FlushHitStack();
1250   }
1251
1252 }
1253
1254 //_____________________________________________________________________________
1255 void AliTRD::RemapTrackHitIDs(Int_t *map)
1256 {
1257   //
1258   // Remap the track IDs
1259   //
1260
1261   if (!fTrackHits) {
1262     return;
1263   }
1264
1265   if (fTrackHits) {
1266     TClonesArray *arr = fTrackHits->GetArray();;
1267     for (Int_t i = 0; i < arr->GetEntriesFast(); i++){
1268       AliTrackHitsParamV2 *info = (AliTrackHitsParamV2 *) (arr->At(i));
1269       info->fTrackID = map[info->fTrackID];
1270     }
1271   }
1272
1273 }
1274
1275 //_____________________________________________________________________________
1276 void AliTRD::ResetHits()
1277 {
1278   //
1279   // Reset the hits
1280   //
1281
1282   AliDetector::ResetHits();
1283   if (fTrackHits) {
1284     fTrackHits->Clear();
1285   }
1286
1287 }
1288
1289 //_____________________________________________________________________________
1290 AliHit* AliTRD::FirstHit(Int_t track)
1291 {
1292   //
1293   // Return the first hit of a track
1294   //
1295
1296   if (fHitType > 1) {
1297     return FirstHit2(track);
1298   }
1299
1300   return AliDetector::FirstHit(track);
1301
1302 }
1303
1304 //_____________________________________________________________________________
1305 AliHit* AliTRD::NextHit()
1306 {
1307   //
1308   // Returns the next hit of a track
1309   //
1310
1311   if (fHitType > 1) {
1312     return NextHit2();
1313   }
1314
1315   return AliDetector::NextHit();
1316
1317 }
1318
1319 //_____________________________________________________________________________
1320 AliHit* AliTRD::FirstHit2(Int_t track) 
1321 {
1322   //
1323   // Initializes the hit iterator.
1324   // Returns the address of the first hit of a track.
1325   // If <track> >= 0 the track is read from disk,
1326   // while if <track> < 0 the first hit of the current
1327   // track is returned.
1328   //
1329
1330   if (track >= 0) {
1331     gAlice->ResetHits();
1332     gAlice->TreeH()->GetEvent(track);
1333   }
1334   
1335   if (fTrackHits) {
1336     fTrackHits->First();
1337     return (AliHit*) fTrackHits->GetHit();
1338   }
1339   else {
1340     return 0;
1341   }
1342
1343 }
1344
1345 //_____________________________________________________________________________
1346 AliHit* AliTRD::NextHit2()
1347 {
1348   //
1349   // Returns the next hit of the current track
1350   //
1351
1352   if (fTrackHits) {
1353     fTrackHits->Next();
1354     return (AliHit *) fTrackHits->GetHit();
1355   }
1356   else {
1357     return 0;
1358   }
1359
1360 }
1361
1362 //_____________________________________________________________________________
1363 void AliTRD::MakeBranch2(Option_t *option, const char *file)
1364 {
1365   //
1366   // Create a new branch in the current Root tree.
1367   // The branch of fHits is automatically split.
1368   //
1369
1370   if (fHitType < 2) {
1371     return;
1372   }
1373
1374   char branchname[10];
1375   sprintf(branchname,"%s2",GetName());
1376
1377   // Get the pointer to the header
1378   const char *cH = strstr(option,"H");
1379  
1380   if (!fTrackHits) {
1381     fTrackHits = new AliTRDtrackHits();
1382   }
1383
1384   if (fTrackHits && gAlice->TreeH() && cH) {
1385
1386     gAlice->TreeH()->Branch(branchname,"AliTRDtrackHits"
1387                                       ,&fTrackHits
1388                                       ,fBufferSize,99);
1389
1390     if (GetDebug() > 1) {
1391       printf("<AliTRD::MakeBranch2> Making Branch %s for trackhits\n"
1392             ,branchname);
1393     }
1394
1395     const char kFolder[] = "RunMC/Event/Data";
1396
1397     if (GetDebug()) {
1398       printf("<AliTRD::MakeBranch2> %15s: Publishing %s to %s\n"
1399             ,ClassName(),branchname,kFolder);
1400     }
1401
1402     Publish(kFolder,&fTrackHits,branchname);
1403
1404     if (file) {
1405       TBranch *b = gAlice->TreeH()->GetBranch(branchname);
1406       TDirectory *wd = gDirectory;
1407       b->SetFile(file);
1408       TIter next(b->GetListOfBranches());
1409       while ((b = (TBranch*) next())) {
1410         b->SetFile(file);
1411       }
1412       wd->cd();
1413       if (GetDebug() > 1) {
1414         printf("<AliTRD::MakeBranch2> Diverting branch %s to file %s\n"
1415               ,branchname,file);
1416       }
1417     }
1418
1419   }
1420
1421 }
1422
1423 //_____________________________________________________________________________
1424 void AliTRD::SetTreeAddress2()
1425 {
1426   //
1427   // Set the branch address for the trackHits tree
1428   //
1429
1430   TBranch *branch;
1431
1432   char branchname[20];
1433
1434   sprintf(branchname,"%s2",GetName());
1435   
1436   // Branch address for hit tree
1437   TTree *treeH = gAlice->TreeH();
1438   if ((treeH) && (fHitType > 0)) {
1439     branch = treeH->GetBranch(branchname);
1440     if (branch) {
1441       branch->SetAddress(&fTrackHits);
1442     }
1443   }
1444
1445 }
1446
1447 //_____________________________________________________________________________
1448 void AliTRD::AddHit2(Int_t track, Int_t det, Float_t *hits, Int_t q
1449                    , Bool_t inDrift)
1450 {
1451   //
1452   // Add a hit to the list
1453   //
1454
1455   Int_t rtrack;
1456
1457   if (fIshunt) {
1458     Int_t primary = gAlice->GetPrimary(track);
1459     gAlice->Particle(primary)->SetBit(kKeepBit);
1460     rtrack = primary;
1461   } 
1462   else {
1463     rtrack = track;
1464     gAlice->FlagTrack(track);
1465   }
1466
1467   if ((fTrackHits) && (fHitType > 0)) {
1468     fTrackHits->AddHitTRD(det,rtrack,hits[0],hits[1],hits[2],q,inDrift);
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
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810