]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/AliTRD.cxx
Changed timebin 0 to be the one closest to the readout
[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.21  2000/06/09 11:10:07  cblume
19 Compiler warnings and coding conventions, next round
20
21 Revision 1.20  2000/06/08 18:32:57  cblume
22 Make code compliant to coding conventions
23
24 Revision 1.19  2000/06/07 16:25:37  cblume
25 Try to remove compiler warnings on Sun and HP
26
27 Revision 1.18  2000/05/08 16:17:27  cblume
28 Merge TRD-develop
29
30 Revision 1.17.2.1  2000/05/08 14:28:59  cblume
31 Introduced SetPHOShole() and SetRICHhole(). AliTRDrecPoint container is now a TObjArray
32
33 Revision 1.17  2000/02/28 19:10:26  cblume
34 Include the new TRD classes
35
36 Revision 1.16.2.2  2000/02/28 17:53:24  cblume
37 Introduce TRD geometry classes
38
39 Revision 1.16.2.1  2000/02/28 17:04:19  cblume
40 Include functions and data members for AliTRDrecPoint
41
42 Revision 1.16  2000/01/19 17:17:35  fca
43 Introducing a list of lists of hits -- more hits allowed for detector now
44
45 Revision 1.15  1999/11/02 17:04:25  fca
46 Small syntax change for HP compiler
47
48 Revision 1.14  1999/11/02 16:57:02  fca
49 Avoid non ansi warnings on HP compilers
50
51 Revision 1.13  1999/11/02 16:35:56  fca
52 New version of TRD introduced
53
54 Revision 1.12  1999/11/01 20:41:51  fca
55 Added protections against using the wrong version of FRAME
56
57 Revision 1.11  1999/09/29 09:24:34  fca
58 Introduction of the Copyright and cvs Log
59
60 */
61
62 ///////////////////////////////////////////////////////////////////////////////
63 //                                                                           //
64 //  Transition Radiation Detector                                            //
65 //  This class contains the basic functions for the Transition Radiation     //
66 //  Detector.                                                                //
67 //                                                                           //
68 ///////////////////////////////////////////////////////////////////////////////
69
70 #include <stdlib.h>
71
72 #include <TMath.h>
73 #include <TNode.h>
74 #include <TPGON.h> 
75 #include <TGeometry.h>
76 #include <TTree.h>
77
78 #include "AliTRD.h"
79 #include "AliRun.h"
80 #include "AliConst.h"
81 #include "AliTRDdigitizer.h"
82 #include "AliTRDclusterizer.h"
83 #include "AliTRDgeometryHole.h"
84 #include "AliTRDgeometryFull.h"
85 #include "AliTRDrecPoint.h"
86 #include "AliMagF.h"
87 #include "AliMC.h"
88  
89 ClassImp(AliTRD)
90  
91 //_____________________________________________________________________________
92 AliTRD::AliTRD()
93 {
94   //
95   // Default constructor
96   //
97
98   fIshunt      = 0;
99   fGasMix      = 0;
100   fHits        = 0;
101   fDigits      = 0;
102
103   fRecPoints   = 0;
104   fNRecPoints  = 0;
105
106   fGeometry    = 0;
107
108 }
109  
110 //_____________________________________________________________________________
111 AliTRD::AliTRD(const char *name, const char *title)
112        : AliDetector(name,title)
113 {
114   //
115   // Standard constructor for the TRD
116   //
117
118   // Check that FRAME is there otherwise we have no place where to
119   // put TRD
120   AliModule* frame = gAlice->GetModule("FRAME");
121   if (!frame) {
122     Error("Ctor","TRD needs FRAME to be present\n");
123     exit(1);
124   } 
125
126   // Define the TRD geometry according to the FRAME geometry
127   if      (frame->IsVersion() == 0) {
128     // Geometry with hole
129     fGeometry = new AliTRDgeometryHole();
130   }
131   else if (frame->IsVersion() == 1) {
132     // Geometry without hole
133     fGeometry = new AliTRDgeometryFull();
134   }
135   else {
136     Error("Ctor","Could not find valid FRAME version\n");
137     exit(1);
138   }
139
140   // Allocate the hit array
141   fHits       = new TClonesArray("AliTRDhit"     ,405);
142   gAlice->AddHitList(fHits);
143
144   // Allocate the digits array
145   fDigits     = 0;
146
147   // Allocate the rec point array
148   fRecPoints  = new TObjArray(400);
149   fNRecPoints = 0;
150    
151   fIshunt = 0;
152   fGasMix = 0;
153
154   SetMarkerColor(kWhite);   
155
156 }
157
158 //_____________________________________________________________________________
159 AliTRD::AliTRD(const AliTRD &trd)
160 {
161   //
162   // Copy constructor
163   //
164
165   ((AliTRD &) trd).Copy(*this);
166
167 }
168
169 //_____________________________________________________________________________
170 AliTRD::~AliTRD()
171 {
172   //
173   // TRD destructor
174   //
175
176   fIshunt = 0;
177
178   delete fGeometry;
179   delete fHits;
180   delete fRecPoints;
181
182 }
183
184 //_____________________________________________________________________________
185 void AliTRD::AddRecPoint(Float_t *pos, Int_t *digits, Int_t det, Float_t amp)
186 {
187   //
188   // Add a reconstructed point for the TRD
189   //
190   
191   AliTRDrecPoint *recPoint = new AliTRDrecPoint();
192   TVector3        posVec(pos[0],pos[1],pos[2]);
193   recPoint->SetLocalPosition(posVec);
194   recPoint->SetDetector(det);
195   recPoint->SetEnergy(amp);
196   for (Int_t iDigit = 0; iDigit < 3; iDigit++) {
197     recPoint->AddDigit(digits[iDigit]);
198   }
199
200   fRecPoints->Add(recPoint);
201
202 }
203
204 //_____________________________________________________________________________
205 void AliTRD::AddDigit(Int_t *digits, Int_t *amp)
206 {
207   //
208   // Add a digit for the TRD
209   //
210
211   TClonesArray &ldigits = *fDigits;
212   new(ldigits[fNdigits++]) AliTRDdigit(kFALSE,digits,amp);
213
214 }
215
216 //_____________________________________________________________________________
217 void AliTRD::AddHit(Int_t track, Int_t *det, Float_t *hits)
218 {
219   //
220   // Add a hit for the TRD
221   //
222
223   TClonesArray &lhits = *fHits;
224   new(lhits[fNhits++]) AliTRDhit(fIshunt,track,det,hits);
225
226 }
227
228 //_____________________________________________________________________________
229 void AliTRD::BuildGeometry()
230 {
231   //
232   // Create the ROOT TNode geometry for the TRD
233   //
234
235   TNode *node, *top;
236   TPGON *pgon;
237   const Int_t kColorTRD = 46;
238   
239   // Find the top node alice
240   top = gAlice->GetGeometry()->GetNode("alice");
241   
242   pgon = new TPGON("S_TRD","TRD","void",0,360,kNsect,4);
243   Float_t ff    = TMath::Cos(kDegrad * 180 / kNsect);
244   Float_t rrmin = kRmin / ff;
245   Float_t rrmax = kRmax / ff;
246   pgon->DefineSection(0,-kZmax1,rrmax,rrmax);
247   pgon->DefineSection(1,-kZmax2,rrmin,rrmax);
248   pgon->DefineSection(2, kZmax2,rrmin,rrmax);
249   pgon->DefineSection(3, kZmax1,rrmax,rrmax);
250   top->cd();
251   node = new TNode("TRD","TRD","S_TRD",0,0,0,"");
252   node->SetLineColor(kColorTRD);
253   fNodes->Add(node);
254
255 }
256  
257 //_____________________________________________________________________________
258 void AliTRD::Copy(TObject &trd)
259 {
260   //
261   // Copy function
262   //
263
264   ((AliTRD &) trd).fGasMix     = fGasMix;
265   ((AliTRD &) trd).fGeometry   = fGeometry;
266   ((AliTRD &) trd).fRecPoints  = fRecPoints;
267   ((AliTRD &) trd).fNRecPoints = fNRecPoints;
268
269   //AliDetector::Copy(trd);
270
271 }
272
273 //_____________________________________________________________________________
274 void AliTRD::CreateGeometry()
275 {
276   //
277   // Creates the volumes for the TRD chambers
278   //
279
280   // Check that FRAME is there otherwise we have no place where to put the TRD
281   AliModule* frame = gAlice->GetModule("FRAME");
282   if (!frame) {
283     printf(" The TRD needs the FRAME to be defined first\n");
284     return;
285   }
286
287   fGeometry->CreateGeometry(fIdtmed->GetArray() - 1299);
288
289 }
290  
291 //_____________________________________________________________________________
292 void AliTRD::CreateMaterials()
293 {
294   //
295   // Create the materials for the TRD
296   // Origin Y.Foka
297   //
298
299   Int_t   isxfld = gAlice->Field()->Integ();
300   Float_t sxmgmx = gAlice->Field()->Max();
301   
302   // For polyethilene (CH2) 
303   Float_t ape[2] = { 12., 1. };
304   Float_t zpe[2] = {  6., 1. };
305   Float_t wpe[2] = {  1., 2. };
306   Float_t dpe    = 0.95;
307
308   // For mylar (C5H4O2) 
309   Float_t amy[3] = { 12., 1., 16. };
310   Float_t zmy[3] = {  6., 1.,  8. };
311   Float_t wmy[3] = {  5., 4.,  2. };
312   Float_t dmy    = 1.39;
313
314   // For CO2 
315   Float_t aco[2] = { 12., 16. };
316   Float_t zco[2] = {  6.,  8. };
317   Float_t wco[2] = {  1.,  2. };
318   Float_t dco    = 0.001977;
319
320   // For water
321   Float_t awa[2] = {  1., 16. };
322   Float_t zwa[2] = {  1.,  8. };
323   Float_t wwa[2] = {  2.,  1. };
324   Float_t dwa    = 1.0;
325
326   // For isobutane (C4H10)
327   Float_t ais[2] = { 12.,  1. };
328   Float_t zis[2] = {  6.,  1. };
329   Float_t wis[2] = {  4., 10. };
330   Float_t dis    = 0.00267;
331
332   // For Xe/CO2-gas-mixture 
333   // Xe-content of the Xe/CO2-mixture (90% / 10%) 
334   Float_t fxc    = .90;
335   // Xe-content of the Xe/Isobutane-mixture (97% / 3%) 
336   Float_t fxi    = .97;
337   Float_t dxe    = .005858;
338   
339   // General tracking parameter
340   Float_t tmaxfd = -10.;
341   Float_t stemax = -1e10;
342   Float_t deemax = -0.1;
343   Float_t epsil  =  1e-4;
344   Float_t stmin  = -0.001;
345   
346   Float_t absl, radl, d, buf[1];
347   Float_t agm[2], dgm, zgm[2], wgm[2];
348   Int_t   nbuf;
349   
350   //////////////////////////////////////////////////////////////////////////
351   //     Define Materials 
352   //////////////////////////////////////////////////////////////////////////
353
354   AliMaterial( 1, "Al $",  26.98, 13.0, 2.7     ,     8.9 ,    37.2);
355   AliMaterial( 2, "Air$",  14.61,  7.3, 0.001205, 30420.0 , 67500.0);
356   AliMaterial( 4, "Xe $", 131.29, 54.0, dxe     ,  1447.59,     0.0);
357   AliMaterial( 5, "Cu $",  63.54, 29.0, 8.96    ,     1.43,    14.8);
358   AliMaterial( 6, "C  $",  12.01,  6.0, 2.265   ,    18.8 ,    74.4);
359   AliMaterial(12, "G10$",  20.00, 10.0, 1.7     ,    19.4 ,   999.0);
360
361   // Mixtures 
362   AliMixture(3, "Polyethilene$",   ape, zpe, dpe, -2, wpe);
363   AliMixture(7, "Mylar$",          amy, zmy, dmy, -3, wmy);
364   AliMixture(8, "CO2$",            aco, zco, dco, -2, wco);
365   AliMixture(9, "Isobutane$",      ais, zis, dis, -2, wis);
366   AliMixture(13,"Water$",          awa, zwa, dwa, -2, wwa);
367
368   // Gas mixtures
369   Char_t namate[21];
370   // Xe/CO2-mixture
371   // Get properties of Xe 
372   gMC->Gfmate((*fIdmate)[4], namate, agm[0], zgm[0], d, radl, absl, buf, nbuf);
373   // Get properties of CO2 
374   gMC->Gfmate((*fIdmate)[8], namate, agm[1], zgm[1], d, radl, absl, buf, nbuf);
375   // Create gas mixture 
376   wgm[0] = fxc;
377   wgm[1] = 1. - fxc;
378   dgm    = wgm[0] * dxe + wgm[1] * dco;
379   AliMixture(10, "Gas mixture 1$", agm, zgm, dgm,  2, wgm);
380   // Xe/Isobutane-mixture
381   // Get properties of Xe 
382   gMC->Gfmate((*fIdmate)[4], namate, agm[0], zgm[0], d, radl, absl, buf, nbuf);
383   // Get properties of Isobutane
384   gMC->Gfmate((*fIdmate)[9], namate, agm[1], zgm[1], d, radl, absl, buf, nbuf);
385   // Create gas mixture 
386   wgm[0] = fxi;
387   wgm[1] = 1. - fxi;
388   dgm    = wgm[0] * dxe + wgm[1] * dis;
389   AliMixture(11, "Gas mixture 2$", agm, zgm, dgm,  2, wgm);
390  
391   //////////////////////////////////////////////////////////////////////////
392   //     Tracking Media Parameters 
393   //////////////////////////////////////////////////////////////////////////
394
395   // Al Frame 
396   AliMedium(1, "Al Frame$",   1, 0, isxfld, sxmgmx
397                 , tmaxfd, stemax, deemax, epsil, stmin);
398   // Air 
399   AliMedium(2, "Air$",        2, 0, isxfld, sxmgmx
400                 , tmaxfd, stemax, deemax, epsil, stmin);
401   // Polyethilene 
402   AliMedium(3, "Radiator$",   3, 0, isxfld, sxmgmx
403                 , tmaxfd, stemax, deemax, epsil, stmin);
404   // Xe 
405   AliMedium(4, "Xe$",         4, 1, isxfld, sxmgmx
406                 , tmaxfd, stemax, deemax, epsil, stmin);
407   // Cu pads 
408   AliMedium(5, "Padplane$",   5, 1, isxfld, sxmgmx
409                 , tmaxfd, stemax, deemax, epsil, stmin);
410   // Fee + cables 
411   AliMedium(6, "Readout$",    1, 0, isxfld, sxmgmx
412                 , tmaxfd, stemax, deemax, epsil, stmin);
413   // C frame 
414   AliMedium(7, "C Frame$",    6, 0, isxfld, sxmgmx
415                 , tmaxfd, stemax, deemax, epsil, stmin);
416   // Mylar foils 
417   AliMedium(8, "Mylar$",      7, 0, isxfld, sxmgmx
418                 , tmaxfd, stemax, deemax, epsil, stmin);
419   if (fGasMix == 1) {
420     // Gas-mixture (Xe/CO2) 
421     AliMedium(9, "Gas-mix$",   10, 1, isxfld, sxmgmx
422                   , tmaxfd, stemax, deemax, epsil, stmin);
423   }
424   else {
425     // Gas-mixture (Xe/Isobutane) 
426     AliMedium(9, "Gas-mix$",   11, 1, isxfld, sxmgmx
427                   , tmaxfd, stemax, deemax, epsil, stmin);
428   }
429   // Nomex-honeycomb (use carbon for the time being) 
430   AliMedium(10, "Nomex$",      6, 0, isxfld, sxmgmx
431                 , tmaxfd, stemax, deemax, epsil, stmin);
432   // Kapton foils (use Mylar for the time being) 
433   AliMedium(11, "Kapton$",     7, 0, isxfld, sxmgmx
434                 , tmaxfd, stemax, deemax, epsil, stmin);
435   // Gas-filling of the radiator 
436   AliMedium(12, "CO2$",        8, 0, isxfld, sxmgmx
437                 , tmaxfd, stemax, deemax, epsil, stmin);
438   // G10-plates
439   AliMedium(13, "G10-plates$",12, 0, isxfld, sxmgmx
440                 , tmaxfd, stemax, deemax, epsil, stmin);
441   // Cooling water
442   AliMedium(14, "Water$",     13, 0, isxfld, sxmgmx
443                 , tmaxfd, stemax, deemax, epsil, stmin);
444
445 }
446
447 //_____________________________________________________________________________
448 void AliTRD::DrawModule()
449 {
450   //
451   // Draw a shaded view of the Transition Radiation Detector version 0
452   //
453
454   // Set everything unseen
455   gMC->Gsatt("*"   ,"SEEN",-1);
456   
457   // Set ALIC mother transparent
458   gMC->Gsatt("ALIC","SEEN", 0);
459   
460   // Set the volumes visible
461   if (fGeometry->IsVersion() == 0) {
462     gMC->Gsatt("B071","SEEN", 0);
463     gMC->Gsatt("B074","SEEN", 0);
464     gMC->Gsatt("B075","SEEN", 0);
465     gMC->Gsatt("B077","SEEN", 0);
466     gMC->Gsatt("BTR1","SEEN", 0);
467     gMC->Gsatt("BTR2","SEEN", 0);
468     gMC->Gsatt("BTR3","SEEN", 0);
469     gMC->Gsatt("TRD1","SEEN", 0);
470     gMC->Gsatt("TRD2","SEEN", 0);
471     gMC->Gsatt("TRD3","SEEN", 0);
472   }
473   else {
474     gMC->Gsatt("B071","SEEN", 0);
475     gMC->Gsatt("B074","SEEN", 0);
476     gMC->Gsatt("B075","SEEN", 0);
477     gMC->Gsatt("B077","SEEN", 0);
478     gMC->Gsatt("BTR1","SEEN", 0);
479     gMC->Gsatt("BTR2","SEEN", 0);
480     gMC->Gsatt("BTR3","SEEN", 0);
481     gMC->Gsatt("TRD1","SEEN", 0);
482     if (fGeometry->GetPHOShole())
483       gMC->Gsatt("TRD2","SEEN", 0);
484     if (fGeometry->GetRICHhole())
485       gMC->Gsatt("TRD3","SEEN", 0);
486   }
487   gMC->Gsatt("UCII","SEEN", 0);
488   gMC->Gsatt("UCIM","SEEN", 0);
489   gMC->Gsatt("UCIO","SEEN", 0);
490   gMC->Gsatt("UL02","SEEN", 1);
491   gMC->Gsatt("UL05","SEEN", 1);
492   gMC->Gsatt("UL06","SEEN", 1);
493   
494   gMC->Gdopt("hide", "on");
495   gMC->Gdopt("shad", "on");
496   gMC->Gsatt("*", "fill", 7);
497   gMC->SetClipBox(".");
498   gMC->SetClipBox("*", 0, 2000, -2000, 2000, -2000, 2000);
499   gMC->DefaultRange();
500   gMC->Gdraw("alic", 40, 30, 0, 12, 9.4, .021, .021);
501   gMC->Gdhead(1111, "Transition Radiation Detector");
502   gMC->Gdman(18, 4, "MAN");
503
504 }
505
506 //_____________________________________________________________________________
507 Int_t AliTRD::DistancetoPrimitive(Int_t , Int_t )
508 {
509   //
510   // Distance between the mouse and the TRD detector on the screen
511   // Dummy routine
512   
513   return 9999;
514
515 }
516  
517 //_____________________________________________________________________________
518 void AliTRD::Init()
519 {
520   //
521   // Initialize the TRD detector after the geometry has been created
522   //
523
524   Int_t i;
525
526   printf("\n");
527   for (i = 0; i < 35; i++) printf("*");
528   printf(" TRD_INIT ");
529   for (i = 0; i < 35; i++) printf("*");
530   printf("\n");
531   printf("\n");
532
533   if      (fGeometry->IsVersion() == 0) {
534     printf("          Geometry for spaceframe with holes initialized.\n\n");
535   }
536   else if (fGeometry->IsVersion() == 1) {
537     printf("          Geometry for spaceframe without holes initialized.\n");
538     if (fGeometry->GetPHOShole())
539       printf("          Leave space in front of PHOS free.\n");
540     if (fGeometry->GetRICHhole())
541       printf("          Leave space in front of RICH free.\n");
542     printf("\n");
543   }
544
545   if (fGasMix == 1)
546     printf("          Gas Mixture: 90%% Xe + 10%% CO2\n\n");
547   else
548     printf("          Gas Mixture: 97%% Xe + 3%% Isobutane\n\n");
549
550 }
551
552 //_____________________________________________________________________________
553 void AliTRD::MakeBranch(Option_t* option)
554 {
555   //
556   // Create Tree branches for the TRD digits and cluster.
557   //
558
559   Int_t  buffersize = 4000;
560   Char_t branchname[15];
561
562   AliDetector::MakeBranch(option);
563
564   Char_t *r = strstr(option,"R");
565   sprintf(branchname,"%srecPoints",GetName());
566   if (fRecPoints && gAlice->TreeR() && r) {
567     gAlice->TreeR()->Branch(branchname,fRecPoints->IsA()->GetName()
568                            ,&fRecPoints,buffersize,0);
569     printf("* AliTRD::MakeBranch * Making Branch %s for points in TreeR\n",branchname);
570   }
571
572 }
573
574 //_____________________________________________________________________________
575 void AliTRD::ResetRecPoints()
576 {
577   //
578   // Reset number of reconstructed points and the point array
579   //
580
581   fNRecPoints = 0;
582   if (fRecPoints) fRecPoints->Delete();
583
584 }
585
586 //_____________________________________________________________________________
587 void AliTRD::SetTreeAddress()
588 {
589   //
590   // Set the branch addresses for the trees.
591   //
592
593   Char_t branchname[15];
594
595   AliDetector::SetTreeAddress();
596
597   TBranch *branch;
598   TTree   *treeR = gAlice->TreeR();
599
600   if (treeR) {
601     sprintf(branchname,"%srecPoints",GetName());
602     if (fRecPoints) {
603       branch = treeR->GetBranch(branchname);
604       if (branch) {
605         branch->SetAddress(&fRecPoints);
606       }
607     }
608   }
609
610 }
611
612 //_____________________________________________________________________________
613 void AliTRD::SetGasMix(Int_t imix)
614 {
615   //
616   // Defines the gas mixture (imix=0:  Xe/Isobutane imix=1: Xe/CO2)
617   //
618   
619   if ((imix < 0) || (imix > 1)) {
620     printf("Wrong input value: %d\n",imix);
621     printf("Use standard setting\n");
622     fGasMix = 0;
623     return;
624   }
625
626   fGasMix = imix;
627
628 }
629
630 //_____________________________________________________________________________
631 AliTRD &AliTRD::operator=(const AliTRD &trd)
632 {
633   //
634   // Assignment operator
635   //
636
637   if (this != &trd) ((AliTRD &) trd).Copy(*this);
638   return *this;
639
640 }