1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.42 2002/10/22 15:53:08 alibrary
19 Introducing Riostream.h
21 Revision 1.41 2002/10/14 14:57:43 hristov
22 Merging the VirtualMC branch to the main development branch (HEAD)
24 Revision 1.36.6.2 2002/07/24 10:09:30 alibrary
27 Revision 1.40 2002/06/13 08:11:56 cblume
28 Add the track references
30 Revision 1.39 2002/06/12 09:54:35 cblume
31 Update of tracking code provided by Sergei
33 Revision 1.38 2002/03/28 14:59:07 cblume
36 Revision 1.37 2002/03/25 20:01:49 cblume
37 Introduce parameter class
39 Revision 1.36 2002/02/11 14:25:27 cblume
40 Geometry update, compressed hit structure
42 Revision 1.35 2001/11/14 12:08:44 cblume
43 Remove unneccessary header files
45 Revision 1.34 2001/11/14 10:50:45 cblume
46 Changes in digits IO. Add merging of summable digits
48 Revision 1.33 2001/11/06 17:19:41 cblume
49 Add detailed geometry and simple simulator
51 Revision 1.32 2001/10/08 06:57:33 hristov
52 Branches for TRD digits are created only during the digitisation
54 Revision 1.31 2001/08/30 09:30:30 hristov
55 The split level of branches is set to 99
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)
60 Revision 1.29 2001/05/21 16:45:47 hristov
61 Last minute changes (C.Blume)
63 Revision 1.28 2001/05/16 14:57:27 alibrary
64 New files for folders and Stack
66 Revision 1.27 2001/05/08 07:05:02 hristov
67 Loop variable declared once (HP, Sun)
69 Revision 1.26 2001/05/07 08:03:22 cblume
70 Generate also hits in the amplification region
72 Revision 1.25 2001/03/13 09:30:35 cblume
73 Update of digitization. Moved digit branch definition to AliTRD
75 Revision 1.24 2001/01/26 19:56:49 hristov
76 Major upgrade of AliRoot code
78 Revision 1.23 2000/11/01 14:53:20 cblume
79 Merge with TRD-develop
81 Revision 1.17.2.6 2000/10/15 23:29:08 cblume
82 Introduced more detailed geometry for the display
84 Revision 1.17.2.5 2000/10/06 16:49:46 cblume
87 Revision 1.17.2.4 2000/10/04 16:34:57 cblume
88 Replace include files by forward declarations
90 Revision 1.17.2.3 2000/09/22 14:45:17 cblume
91 Included changes for the tracking
93 Revision 1.17.2.2 2000/09/18 13:25:13 cblume
94 Included LoadPoints() method to display the TR photons
96 Revision 1.22 2000/10/02 21:28:19 fca
97 Removal of useless dependecies via forward declarations
99 Revision 1.21 2000/06/09 11:10:07 cblume
100 Compiler warnings and coding conventions, next round
102 Revision 1.20 2000/06/08 18:32:57 cblume
103 Make code compliant to coding conventions
105 Revision 1.19 2000/06/07 16:25:37 cblume
106 Try to remove compiler warnings on Sun and HP
108 Revision 1.18 2000/05/08 16:17:27 cblume
111 Revision 1.21 2000/06/09 11:10:07 cblume
112 Compiler warnings and coding conventions, next round
114 Revision 1.20 2000/06/08 18:32:57 cblume
115 Make code compliant to coding conventions
117 Revision 1.19 2000/06/07 16:25:37 cblume
118 Try to remove compiler warnings on Sun and HP
120 Revision 1.18 2000/05/08 16:17:27 cblume
123 Revision 1.17.2.1 2000/05/08 14:28:59 cblume
124 Introduced SetPHOShole() and SetRICHhole(). AliTRDrecPoint container is now a TObjArray
126 Revision 1.17 2000/02/28 19:10:26 cblume
127 Include the new TRD classes
129 Revision 1.16.2.2 2000/02/28 17:53:24 cblume
130 Introduce TRD geometry classes
132 Revision 1.16.2.1 2000/02/28 17:04:19 cblume
133 Include functions and data members for AliTRDrecPoint
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
138 Revision 1.15 1999/11/02 17:04:25 fca
139 Small syntax change for HP compiler
141 Revision 1.14 1999/11/02 16:57:02 fca
142 Avoid non ansi warnings on HP compilers
144 Revision 1.13 1999/11/02 16:35:56 fca
145 New version of TRD introduced
147 Revision 1.12 1999/11/01 20:41:51 fca
148 Added protections against using the wrong version of FRAME
150 Revision 1.11 1999/09/29 09:24:34 fca
151 Introduction of the Copyright and cvs Log
155 ///////////////////////////////////////////////////////////////////////////////
157 // Transition Radiation Detector //
158 // This class contains the basic functions for the Transition Radiation //
161 ///////////////////////////////////////////////////////////////////////////////
164 #include <Riostream.h>
168 #include <TGeometry.h>
173 #include <TParticle.h>
174 #include <TLorentzVector.h>
177 #include "AliConst.h"
178 #include "AliDigit.h"
181 #include "AliTrackReference.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"
198 //_____________________________________________________________________________
202 // Default constructor
226 //_____________________________________________________________________________
227 AliTRD::AliTRD(const char *name, const char *title)
228 : AliDetector(name,title)
231 // Standard constructor for the TRD
234 // Check that FRAME is there otherwise we have no place where to
236 AliModule* frame = gAlice->GetModule("FRAME");
238 Error("Ctor","TRD needs FRAME to be present\n");
242 // Define the TRD geometry according to the FRAME geometry
243 if (frame->IsVersion() == 0) {
244 // Geometry with hole
245 fGeometry = new AliTRDgeometryHole();
247 else if (frame->IsVersion() == 1) {
248 // Geometry without hole
249 fGeometry = new AliTRDgeometryFull();
252 Error("Ctor","Could not find valid FRAME version\n");
256 // Allocate the hit array
257 fHits = new TClonesArray("AliTRDhit" ,405);
258 gAlice->AddHitList(fHits);
260 // Allocate the digits array
263 // Allocate the rec point array
264 fRecPoints = new TObjArray(400);
279 SetMarkerColor(kWhite);
283 //_____________________________________________________________________________
284 AliTRD::AliTRD(const AliTRD &trd)
290 ((AliTRD &) trd).Copy(*this);
294 //_____________________________________________________________________________
322 //_____________________________________________________________________________
323 void AliTRD::AddCluster(Float_t *pos, Int_t det, Float_t amp
324 , Int_t *tracks, Float_t *sig, Int_t iType)
327 // Add a cluster for the TRD
330 AliTRDcluster *c = new AliTRDcluster();
333 c->AddTrackIndex(tracks);
337 c->SetSigmaY2(sig[0]);
338 c->SetSigmaZ2(sig[1]);
339 c->SetLocalTimeBin(((Int_t) pos[2]));
363 //_____________________________________________________________________________
364 void AliTRD::AddTrackReference(Int_t label, TLorentzVector p, TLorentzVector x)
367 // Add a trackrefernce to the list
370 if (!fTrackReferences) {
371 Error("AddTrackReference","Container fTrackRefernce not active\n");
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);
384 //_____________________________________________________________________________
385 void AliTRD::Hits2Digits()
391 AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
392 ,"TRD digitizer class");
393 digitizer->SetDebug(GetDebug());
394 digitizer->SetEvent(gAlice->GetEvNumber());
397 digitizer->InitDetector();
400 digitizer->MakeDigits();
402 // Write the digits into the input file
403 if (digitizer->MakeBranch(fDigitsFile)) {
405 digitizer->WriteDigits();
407 // Save the digitizer class in the AliROOT
414 //_____________________________________________________________________________
415 void AliTRD::Hits2SDigits()
418 // Create summable digits
421 AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
422 ,"TRD digitizer class");
423 digitizer->SetDebug(GetDebug());
425 // For the summable digits
426 digitizer->SetSDigits(kTRUE);
427 digitizer->SetEvent(gAlice->GetEvNumber());
430 digitizer->InitDetector();
432 // Create the TRD s-digits branch
433 digitizer->MakeDigits();
435 // Write the digits into the input file
436 if (digitizer->MakeBranch(fDigitsFile)) {
438 digitizer->WriteDigits();
440 // Save the digitizer class in the AliROOT
447 //_____________________________________________________________________________
448 void AliTRD::SDigits2Digits()
451 // Create final digits from summable digits
454 // Create the TRD digitizer
455 AliTRDdigitizer *digitizer = new AliTRDdigitizer("TRDdigitizer"
456 ,"TRD digitizer class");
457 digitizer->SetDebug(GetDebug());
460 digitizer->SetEvent(gAlice->GetEvNumber());
463 digitizer->InitDetector();
465 // Read the s-digits via digits manager
466 AliTRDdigitsManager *sdigitsManager = new AliTRDdigitsManager();
467 sdigitsManager->SetDebug(GetDebug());
468 sdigitsManager->SetSDigits(kTRUE);
470 sdigitsManager->Open(fDigitsFile);
472 sdigitsManager->CreateArrays();
473 sdigitsManager->ReadDigits();
475 // Add the s-digits to the input list
476 digitizer->AddSDigitsManager(sdigitsManager);
478 // Convert the s-digits to normal digits
479 digitizer->SDigits2Digits();
482 if (digitizer->MakeBranch(fDigitsFile)) {
484 digitizer->WriteDigits();
490 //_____________________________________________________________________________
491 void AliTRD::AddHit(Int_t track, Int_t det, Float_t *hits, Int_t q
495 // Add a hit for the TRD
497 // The data structure is set according to fHitType:
498 // bit0: standard TClonesArray
499 // bit1: compressed trackHits structure
503 TClonesArray &lhits = *fHits;
504 new(lhits[fNhits++]) AliTRDhit(fIshunt,track,det,hits,q);
508 AddHit2(track,det,hits,q,inDrift);
513 //_____________________________________________________________________________
514 void AliTRD::BuildGeometry()
517 // Create the ROOT TNode geometry for the TRD
524 Float_t zmax1, zmax2;
528 const Int_t kColorTRD = 46;
530 // Find the top node alice
531 top = gAlice->GetGeometry()->GetNode("alice");
533 if (fDisplayType == 0) {
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);
543 node = new TNode("TRD","TRD","S_TRD",0,0,0,"");
544 node->SetLineColor(kColorTRD);
548 else if (fDisplayType == 1) {
552 Float_t slope = (AliTRDgeometry::Zmax1() - AliTRDgeometry::Zmax2())
553 / (AliTRDgeometry::Rmax() - AliTRDgeometry::Rmin());
555 rmin = AliTRDgeometry::Rmin() + AliTRDgeometry::CraHght();
556 rmax = rmin + AliTRDgeometry::CdrHght();
558 Float_t thickness = rmin - AliTRDgeometry::Rmin();
559 zmax2 = AliTRDgeometry::Zmax2() + slope * thickness;
560 zmax1 = zmax2 + slope * AliTRDgeometry::DrThick();
562 for (iPlan = 0; iPlan < AliTRDgeometry::Nplan(); iPlan++) {
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);
571 node = new TNode("TRD","TRD",name,0,0,0,"");
572 node->SetLineColor(kColorTRD);
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;
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();
589 for (iPlan = 0; iPlan < AliTRDgeometry::Nplan(); iPlan++) {
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);
598 node = new TNode("TRD","TRD",name,0,0,0,"");
599 node->SetLineColor(kColorTRD);
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;
614 //_____________________________________________________________________________
615 void AliTRD::Copy(TObject &trd)
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;
631 //AliDetector::Copy(trd);
635 //_____________________________________________________________________________
636 void AliTRD::CreateGeometry()
639 // Creates the volumes for the TRD chambers
642 // Check that FRAME is there otherwise we have no place where to put the TRD
643 AliModule* frame = gAlice->GetModule("FRAME");
645 printf(" The TRD needs the FRAME to be defined first\n");
649 fGeometry->CreateGeometry(fIdtmed->GetArray() - 1299);
653 //_____________________________________________________________________________
654 void AliTRD::CreateMaterials()
657 // Create the materials for the TRD
661 Int_t isxfld = gAlice->Field()->Integ();
662 Float_t sxmgmx = gAlice->Field()->Max();
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. };
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. };
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;
683 Float_t awa[2] = { 1., 16. };
684 Float_t zwa[2] = { 1., 8. };
685 Float_t wwa[2] = { 2., 1. };
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;
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 };
700 // For Xe/CO2-gas-mixture
701 // Xe-content of the Xe/CO2-mixture (85% / 15%)
703 // Xe-content of the Xe/Isobutane-mixture (97% / 3%)
705 Float_t dxe = .005858;
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;
714 Float_t absl, radl, d, buf[1];
715 Float_t agm[2], zgm[2], wgm[2];
719 //////////////////////////////////////////////////////////////////////////
721 //////////////////////////////////////////////////////////////////////////
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);
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);
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
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
761 dgm2 = wgm[0] * dxe + wgm[1] * dis;
762 AliMixture(11, "Gas mixture 2", agm, zgm, dgm2, 2, wgm);
764 //////////////////////////////////////////////////////////////////////////
765 // Tracking Media Parameters
766 //////////////////////////////////////////////////////////////////////////
769 AliMedium(1, "Al Frame", 1, 0, isxfld, sxmgmx
770 , tmaxfd, stemax, deemax, epsil, stmin);
772 AliMedium(2, "Air", 2, 0, isxfld, sxmgmx
773 , tmaxfd, stemax, deemax, epsil, stmin);
775 AliMedium(3, "Radiator", 3, 0, isxfld, sxmgmx
776 , tmaxfd, stemax, deemax, epsil, stmin);
778 AliMedium(4, "Xe", 4, 1, isxfld, sxmgmx
779 , tmaxfd, stemax, deemax, epsil, stmin);
781 AliMedium(5, "Padplane", 5, 1, isxfld, sxmgmx
782 , tmaxfd, stemax, deemax, epsil, stmin);
784 AliMedium(6, "Readout", 1, 0, isxfld, sxmgmx
785 , tmaxfd, stemax, deemax, epsil, stmin);
787 AliMedium(7, "C Frame", 6, 0, isxfld, sxmgmx
788 , tmaxfd, stemax, deemax, epsil, stmin);
790 AliMedium(8, "Mylar", 7, 0, isxfld, sxmgmx
791 , tmaxfd, stemax, deemax, epsil, stmin);
793 // Gas-mixture (Xe/CO2)
794 AliMedium(9, "Gas-mix", 10, 1, isxfld, sxmgmx
795 , tmaxfd, stemax, deemax, epsil, stmin);
798 // Gas-mixture (Xe/Isobutane)
799 AliMedium(9, "Gas-mix", 11, 1, isxfld, sxmgmx
800 , tmaxfd, stemax, deemax, epsil, stmin);
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);
812 AliMedium(13, "G10-plates",12, 0, isxfld, sxmgmx
813 , tmaxfd, stemax, deemax, epsil, stmin);
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);
821 AliMedium(16, "MCM-Al" , 1, 0, isxfld, sxmgmx
822 , tmaxfd, stemax, deemax, epsil, stmin);
824 AliMedium(17, "MCM-Sn" , 15, 0, isxfld, sxmgmx
825 , tmaxfd, stemax, deemax, epsil, stmin);
827 AliMedium(18, "MCM-Cu" , 5, 0, isxfld, sxmgmx
828 , tmaxfd, stemax, deemax, epsil, stmin);
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);
839 AliMedium(22, "Conn-PE" , 3, 0, isxfld, sxmgmx
840 , tmaxfd, stemax, deemax, epsil, stmin);
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);
848 // Save the density values for the TRD absorbtion
857 //_____________________________________________________________________________
858 void AliTRD::DrawModule() const
861 // Draw a shaded view of the Transition Radiation Detector version 0
864 // Set everything unseen
865 gMC->Gsatt("*" ,"SEEN",-1);
867 // Set ALIC mother transparent
868 gMC->Gsatt("ALIC","SEEN", 0);
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);
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);
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);
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);
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");
916 //_____________________________________________________________________________
917 Int_t AliTRD::DistancetoPrimitive(Int_t , Int_t ) const
920 // Distance between the mouse and the TRD detector on the screen
927 //_____________________________________________________________________________
931 // Initialize the TRD detector after the geometry has been created
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("*");
944 if (fGeometry->IsVersion() == 0) {
945 printf("%s: Geometry for spaceframe with holes initialized\n",ClassName());
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());
956 printf("%s: Gas Mixture: 85%% Xe + 15%% CO2\n",ClassName());
959 printf("%s: Gas Mixture: 97%% Xe + 3%% Isobutane\n",ClassName());
964 //_____________________________________________________________________________
965 void AliTRD::LoadPoints(Int_t track)
968 // Store x, y, z of all hits in memory.
969 // Hit originating from TR photons are given a different color
973 // AliDetector::LoadPoints(track);
977 if ((fHits == 0) && (fTrackHits == 0)) return;
981 nhits = fHits->GetEntriesFast();
984 nhits = fTrackHits->GetEntriesFast();
986 if (nhits == 0) return;
988 Int_t tracks = gAlice->GetNtrack();
989 if (fPoints == 0) fPoints = new TObjArray(tracks);
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++) {
1008 AliTRDpoints *points = 0;
1011 Int_t chunk = nhits / 4 + 1;
1013 // Loop over all the hits and store their position
1014 ahit = (AliTRDhit *) FirstHit(-1);
1018 if (ahit->GetCharge() >= 0) {
1020 trk = ahit->GetTrack();
1021 if (ntrkE[trk] == limiE[trk]) {
1022 // Initialise a new track
1023 fp = new Float_t[3*(limiE[trk]+chunk)];
1025 memcpy(fp,coorE[trk],sizeof(Float_t)*3*limiE[trk]);
1026 delete [] coorE[trk];
1028 limiE[trk] += chunk;
1034 fp[3*ntrkE[trk] ] = ahit->X();
1035 fp[3*ntrkE[trk]+1] = ahit->Y();
1036 fp[3*ntrkE[trk]+2] = ahit->Z();
1041 else if ((ahit->GetCharge() < 0) && (fDrawTR)) {
1043 trk = ahit->GetTrack();
1044 if (ntrkT[trk] == limiT[trk]) {
1045 // Initialise a new track
1046 fp = new Float_t[3*(limiT[trk]+chunk)];
1048 memcpy(fp,coorT[trk],sizeof(Float_t)*3*limiT[trk]);
1049 delete [] coorT[trk];
1051 limiT[trk] += chunk;
1057 fp[3*ntrkT[trk] ] = ahit->X();
1058 fp[3*ntrkT[trk]+1] = ahit->Y();
1059 fp[3*ntrkT[trk]+2] = ahit->Z();
1064 ahit = (AliTRDhit *) NextHit();
1068 for (trk = 0; trk < tracks; ++trk) {
1070 if (ntrkE[trk] || ntrkT[trk]) {
1072 points = new AliTRDpoints();
1073 points->SetDetector(this);
1074 points->SetParticle(trk);
1076 // Set the dEdx points
1078 points->SetMarkerColor(GetMarkerColor());
1079 points->SetMarkerSize(GetMarkerSize());
1080 points->SetPolyMarker(ntrkE[trk],coorE[trk],GetMarkerStyle());
1081 delete [] coorE[trk];
1085 // Set the TR photon points
1087 points->SetTRpoints(ntrkT[trk],coorT[trk]);
1088 delete [] coorT[trk];
1092 fPoints->AddAt(points,trk);
1107 //_____________________________________________________________________________
1108 void AliTRD::MakeBranch(Option_t* option, const char *file)
1111 // Create Tree branches for the TRD digits.
1114 Int_t buffersize = 4000;
1115 Char_t branchname[15];
1116 sprintf(branchname,"%s",GetName());
1118 const char *cD = strstr(option,"D");
1120 AliDetector::MakeBranch(option,file);
1122 if (fDigits && gAlice->TreeD() && cD) {
1123 MakeBranchInTree(gAlice->TreeD(),branchname,&fDigits,buffersize,file);
1127 MakeBranch2(option,file);
1132 //_____________________________________________________________________________
1133 void AliTRD::ResetDigits()
1136 // Reset number of digits and the digits array for this detector
1140 if (fDigits) fDigits->Clear();
1144 //_____________________________________________________________________________
1145 void AliTRD::ResetRecPoints()
1148 // Reset number of reconstructed points and the point array
1153 Int_t nentr = fRecPoints->GetEntriesFast();
1154 for (Int_t i = 0; i < nentr; i++) delete fRecPoints->RemoveAt(i);
1159 //_____________________________________________________________________________
1160 void AliTRD::SetTreeAddress()
1163 // Set the branch addresses for the trees.
1166 Char_t branchname[15];
1168 AliDetector::SetTreeAddress();
1171 TTree *treeR = gAlice->TreeR();
1174 sprintf(branchname,"%scluster",GetName());
1176 branch = treeR->GetBranch(branchname);
1178 branch->SetAddress(&fRecPoints);
1189 //_____________________________________________________________________________
1190 void AliTRD::SetGasMix(Int_t imix)
1193 // Defines the gas mixture (imix=0: Xe/Isobutane imix=1: Xe/CO2)
1196 if ((imix < 0) || (imix > 1)) {
1197 printf("Wrong input value: %d\n",imix);
1198 printf("Use standard setting\n");
1207 //_____________________________________________________________________________
1208 void AliTRD::SetPHOShole()
1211 // Selects a geometry with a hole in front of the PHOS
1214 fGeometry->SetPHOShole();
1218 //_____________________________________________________________________________
1219 void AliTRD::SetRICHhole()
1222 // Selects a geometry with a hole in front of the RICH
1225 fGeometry->SetRICHhole();
1229 //_____________________________________________________________________________
1230 AliTRD &AliTRD::operator=(const AliTRD &trd)
1233 // Assignment operator
1236 if (this != &trd) ((AliTRD &) trd).Copy(*this);
1241 //_____________________________________________________________________________
1242 void AliTRD::FinishPrimary()
1245 // Store the hits in the containers after all primaries are finished
1249 fTrackHits->FlushHitStack();
1254 //_____________________________________________________________________________
1255 void AliTRD::RemapTrackHitIDs(Int_t *map)
1258 // Remap the track IDs
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];
1275 //_____________________________________________________________________________
1276 void AliTRD::ResetHits()
1282 AliDetector::ResetHits();
1284 fTrackHits->Clear();
1289 //_____________________________________________________________________________
1290 AliHit* AliTRD::FirstHit(Int_t track)
1293 // Return the first hit of a track
1297 return FirstHit2(track);
1300 return AliDetector::FirstHit(track);
1304 //_____________________________________________________________________________
1305 AliHit* AliTRD::NextHit()
1308 // Returns the next hit of a track
1315 return AliDetector::NextHit();
1319 //_____________________________________________________________________________
1320 AliHit* AliTRD::FirstHit2(Int_t track)
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.
1331 gAlice->ResetHits();
1332 gAlice->TreeH()->GetEvent(track);
1336 fTrackHits->First();
1337 return (AliHit*) fTrackHits->GetHit();
1345 //_____________________________________________________________________________
1346 AliHit* AliTRD::NextHit2()
1349 // Returns the next hit of the current track
1354 return (AliHit *) fTrackHits->GetHit();
1362 //_____________________________________________________________________________
1363 void AliTRD::MakeBranch2(Option_t *option, const char *file)
1366 // Create a new branch in the current Root tree.
1367 // The branch of fHits is automatically split.
1374 char branchname[10];
1375 sprintf(branchname,"%s2",GetName());
1377 // Get the pointer to the header
1378 const char *cH = strstr(option,"H");
1381 fTrackHits = new AliTRDtrackHits();
1384 if (fTrackHits && gAlice->TreeH() && cH) {
1386 gAlice->TreeH()->Branch(branchname,"AliTRDtrackHits"
1390 if (GetDebug() > 1) {
1391 printf("<AliTRD::MakeBranch2> Making Branch %s for trackhits\n"
1395 const char kFolder[] = "RunMC/Event/Data";
1398 printf("<AliTRD::MakeBranch2> %15s: Publishing %s to %s\n"
1399 ,ClassName(),branchname,kFolder);
1402 Publish(kFolder,&fTrackHits,branchname);
1405 TBranch *b = gAlice->TreeH()->GetBranch(branchname);
1406 TDirectory *wd = gDirectory;
1408 TIter next(b->GetListOfBranches());
1409 while ((b = (TBranch*) next())) {
1413 if (GetDebug() > 1) {
1414 printf("<AliTRD::MakeBranch2> Diverting branch %s to file %s\n"
1423 //_____________________________________________________________________________
1424 void AliTRD::SetTreeAddress2()
1427 // Set the branch address for the trackHits tree
1432 char branchname[20];
1434 sprintf(branchname,"%s2",GetName());
1436 // Branch address for hit tree
1437 TTree *treeH = gAlice->TreeH();
1438 if ((treeH) && (fHitType > 0)) {
1439 branch = treeH->GetBranch(branchname);
1441 branch->SetAddress(&fTrackHits);
1447 //_____________________________________________________________________________
1448 void AliTRD::AddHit2(Int_t track, Int_t det, Float_t *hits, Int_t q
1452 // Add a hit to the list
1458 Int_t primary = gAlice->GetPrimary(track);
1459 gAlice->Particle(primary)->SetBit(kKeepBit);
1464 gAlice->FlagTrack(track);
1467 if ((fTrackHits) && (fHitType > 0)) {
1468 fTrackHits->AddHitTRD(det,rtrack,hits[0],hits[1],hits[2],q,inDrift);