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