]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv2.cxx
Minor corrections thanks to I.Hrivnacova
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv2.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.8  1999/09/29 09:24:25  fca
19 Introduction of the Copyright and cvs Log
20
21 */
22
23 //-*-C++-*-
24 //_________________________________________________________________________
25 // Manager and hit classes for PHOS
26 //*-- Author : Maxim Volkov, RRC KI
27 // AliPHOSv2 derives directly from AliDetector, because too much functionality
28 // has been put in AliPHOS for my liking.
29 //////////////////////////////////////////////////////////////////////////////
30
31 // --- ROOT system ---
32 #include "TH1.h"
33 #include "TRandom.h"
34 #include "TFile.h"
35 #include "TTree.h"
36 #include "TBRIK.h"
37 #include "TNode.h"
38
39 // --- Standard library ---
40 #include <stdio.h>
41 #include <string.h>
42 #include <stdlib.h>
43
44 // --- galice header files ---
45 #include "AliPHOSv2.h"
46 #include "AliRun.h"
47 #include "AliConst.h"
48
49 ///////////////////////////////////////////////////////////////////////////////
50
51 ClassImp(AliPHOSv2)
52
53 ///////////////////////////////////////////////////////////////////////////////
54   
55 AliPHOSv2::~AliPHOSv2(void)
56 {
57
58   //  fNV       = 0;
59   //  fNH       = 0;
60   fIshunt   = 0;
61   
62   delete fHits;
63
64 }
65
66 ///////////////////////////////////////////////////////////////////////////////
67
68 AliPHOSv2::AliPHOSv2()
69 {
70
71   //  fNV       = 0;
72   //  fNH       = 0;
73   fIshunt   = 0;
74   fHits     = 0;
75   fDigits   = 0;
76
77 }
78
79 ///////////////////////////////////////////////////////////////////////////////
80
81 AliPHOSv2::AliPHOSv2(const char *name, const char *title):
82   AliPHOS(name,title)
83 {
84   
85   //Begin_Html
86   /*
87     <img src="picts/aliphos.gif">
88   */
89   //End_Html
90   
91   fHits   = new TClonesArray("AliPHOShitv2",  405);
92   
93   //  fNV         =  4;
94   //  fNH         =  4;
95   fIshunt     =  1; // All hits are associated with primary particles
96
97   DefPars(); // Set geometry parameters
98   
99   //  SetMarkerColor(kGreen);
100   //  SetMarkerStyle(2);
101   //  SetMarkerSize(0.4);
102   
103 }
104
105 //////////////////////////////////////////////////////////////////////////////
106
107 void AliPHOSv2::DefPars(void)
108 {
109
110   // Initialization of GEANT3 geometry parameters
111   fXtlSize[0]=2.2;
112   fXtlSize[1]=18.0;
113   fXtlSize[2]=2.2;
114
115   fWrapThickness=0.01;
116
117   fPINSize[0]=1.0;
118   fPINSize[1]=0.1;
119   fPINSize[2]=1.0;
120
121   fCPVThickness=1.0;
122
123   fPHOSFoam[0]=214.6;
124   fPHOSFoam[1]=80.0;
125   fPHOSFoam[2]=260.0;
126
127   fPHOStxwall[0]=209.0;
128   fPHOStxwall[1]=71.0;
129   fPHOStxwall[2]=250.0;
130
131   fPHOSAir[0]=206.0;
132   fPHOSAir[1]=66.0;
133   fPHOSAir[2]=244.0;
134
135   fRadius[0]=447.0;
136   fRadius[1]=460.0;
137
138   fPHOSextra[0]=0.005; // Titanium cover thickness
139   fPHOSextra[1]=6.95; // Crystal support height
140   fPHOSextra[2]=4.0; // Thermo Insulating outer cover Upper plate thickness
141   fPHOSextra[3]=5.0; // Upper Polystyrene Foam plate thickness
142   fPHOSextra[4]=2.0; // Thermo insulating Crystal Block wall thickness
143   fPHOSextra[5]=0.06; // Upper Cooling Plate thickness
144   fPHOSextra[6]=10.0; // Al Support Plate thickness
145   fPHOSextra[7]=3.0; // Lower Thermo Insulating Plate thickness
146   fPHOSextra[8]=1.0; // Lower Textolit Plate thickness
147   fPHOSextra[9]=0.03; // 1/2 total gap between adjacent crystals
148
149   fNphi=88;
150   fNz=104;
151
152   fNModules=4;
153
154   fPHOSAngle[0]=0.0; // Module position angles are set in CreateGeometry()
155   fPHOSAngle[1]=0.0;
156   fPHOSAngle[2]=0.0;
157   fPHOSAngle[3]=0.0;
158
159 }
160
161 //////////////////////////////////////////////////////////////////////////////
162
163 void AliPHOSv2::Init(void)
164 {
165
166   Int_t i;
167
168   printf("\n");
169   for(i=0;i<35;i++) printf("*");
170   printf(" PHOS_INIT ");
171   for(i=0;i<35;i++) printf("*");
172   printf("\n");
173
174   // Here the PHOS initialisation code (if any!)
175   for(i=0;i<80;i++) printf("*");
176   printf("\n");
177
178 }
179
180 //////////////////////////////////////////////////////////////////////////////
181
182 void AliPHOSv2::BuildGeometry()
183 {
184
185   // Stolen completely from A. Zvyagine
186
187   TNode *Node, *Top;
188
189   const int kColorPHOS = kRed;
190
191   Top=gAlice->GetGeometry()->GetNode("alice");
192
193   // PHOS
194   Float_t pphi=12.9399462;
195   new TRotMatrix("rot988","rot988",90,-3*pphi,90,90-3*pphi,0,0);
196   new TRotMatrix("rot989","rot989",90,-  pphi,90,90-  pphi,0,0);
197   new TRotMatrix("rot990","rot990",90,   pphi,90,90+  pphi,0,0);
198   new TRotMatrix("rot991","rot991",90, 3*pphi,90,90+3*pphi,0,0);
199   new TBRIK("S_PHOS","PHOS box","void",107.3,40,130);
200   Top->cd();
201   Node=new TNode("PHOS1","PHOS1","S_PHOS",-317.824921,-395.014343,0,"rot988");
202   Node->SetLineColor(kColorPHOS);
203   fNodes->Add(Node);
204   Top->cd();
205   Node=new TNode("PHOS2","PHOS2","S_PHOS",-113.532333,-494.124908,0,"rot989");
206   fNodes->Add(Node);
207   Node->SetLineColor(kColorPHOS);
208   Top->cd();
209   Node=new TNode("PHOS3","PHOS3","S_PHOS", 113.532333,-494.124908,0,"rot990");
210   Node->SetLineColor(kColorPHOS);
211   fNodes->Add(Node);
212   Top->cd();
213   Node=new TNode("PHOS4","PHOS4","S_PHOS", 317.824921,-395.014343,0,"rot991");
214   Node->SetLineColor(kColorPHOS);
215   fNodes->Add(Node);
216
217 }
218
219 //////////////////////////////////////////////////////////////////////////////
220
221 void AliPHOSv2::CreateMaterials()
222 {
223
224   // DEFINITION OF AVAILABLE PHOS MATERIALS
225   
226   Int_t   ISXFLD=gAlice->Field()->Integ();
227   Float_t SXMGMX=gAlice->Field()->Max();
228
229   // --- The PbWO4 crystals ---
230   Float_t AX[3]={207.19, 183.85, 16.0};
231   Float_t ZX[3]={82.0, 74.0, 8.0};
232   Float_t WX[3]={1.0, 1.0, 4.0};
233   Float_t DX=8.28;
234
235   // --- Titanium ---
236   Float_t ATIT[3]={47.88, 26.98, 54.94};
237   Float_t ZTIT[3]={22.0, 13.0, 25.0};
238   Float_t WTIT[3]={69.0, 6.0, 1.0};
239   Float_t DTIT=4.5;
240
241   // --- The polysterene scintillator (CH) ---
242   Float_t AP[2]={12.011, 1.00794};
243   Float_t ZP[2]={6.0, 1.0};
244   Float_t WP[2]={1.0, 1.0};
245   Float_t DP=1.032;
246
247   // --- Tyvek (CnH2n) ---
248   Float_t AT[2]={12.011, 1.00794};
249   Float_t ZT[2]={6.0, 1.0};
250   Float_t WT[2]={1.0, 2.0};
251   Float_t DT=0.331;
252
253   // --- Polystyrene foam ---
254   Float_t AF[2]={12.011, 1.00794};
255   Float_t ZF[2]={6.0, 1.0};
256   Float_t WF[2]={1.0, 1.0};
257   Float_t DF=0.12;
258
259   // --- Foam thermo insulation ---
260   Float_t ATI[2]={12.011, 1.00794};
261   Float_t ZTI[2]={6.0, 1.0};
262   Float_t WTI[2]={1.0, 1.0};
263   Float_t DTI=0.1;
264
265   // --- Textolit ---
266   Float_t ATX[4]={16.0, 28.09, 12.011, 1.00794};
267   Float_t ZTX[4]={8.0, 14.0, 6.0, 1.0};
268   Float_t WTX[4]={292.0, 68.0, 462.0, 736.0};
269   Float_t DTX=1.75;
270
271   Int_t *idtmed = fIdtmed->GetArray()-699;
272
273   AliMixture(0, "PbWO4$", AX, ZX, DX, -3, WX);
274   AliMixture(1, "Polystyrene$", AP, ZP, DP, -2, WP);
275   AliMaterial(2, "Al$", 26.98, 13., 2.7, 8.9, 999., 0, 0);
276   // ---                     Absorption length^ is ignored ---
277   AliMixture(3, "Tyvek$", AT, ZT, DT, -2, WT);
278   AliMixture(4, "Foam$", AF, ZF, DF, -2, WF);
279   AliMixture(5, "Titanium$", ATIT, ZTIT, DTIT, -3, WTIT);
280   AliMaterial(6, "Si$", 28.09, 14., 2.33, 9.36, 42.3, 0, 0);
281   AliMixture(7, "Thermo Insul.$", ATI, ZTI, DTI, -2, WTI);
282   AliMixture(8, "Textolit$", ATX, ZTX, DTX, -4, WTX);
283   AliMaterial(99, "Air$", 14.61, 7.3, 0.001205, 30420., 67500., 0, 0);
284   
285   AliMedium(0, "PHOS Xtal    $", 0, 1,
286             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0);
287   AliMedium(1, "CPV scint.   $", 1, 1,
288             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0);
289   AliMedium(2, "Al parts     $", 2, 0,
290             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0);
291   AliMedium(3, "Tyvek wrapper$", 3, 0,
292             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.001, 0.001, 0, 0);
293   AliMedium(4, "Polyst. foam $", 4, 0,
294             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0);
295   AliMedium(5, "Titan. cover $", 5, 0,
296             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.0001, 0.0001, 0, 0);
297   AliMedium(6, "Si PIN       $", 6, 0,
298             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.01, 0.01, 0, 0);
299   AliMedium(7, "Thermo Insul.$", 7, 0,
300             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0);
301   AliMedium(8, "Textolit     $", 8, 0,
302             ISXFLD, SXMGMX, 10.0, 0.1, 0.1, 0.1, 0.1, 0, 0);
303   AliMedium(99, "Air          $", 99, 0,
304             ISXFLD, SXMGMX, 10.0, 1.0, 0.1, 0.1, 10.0, 0, 0);
305
306   // --- Set decent energy thresholds for gamma and electron tracking
307   gMC->Gstpar(idtmed[699],"CUTGAM",0.5E-4);
308   gMC->Gstpar(idtmed[699],"CUTELE",1.0E-4);
309   // --- Generate explicitly delta rays in the titan cover ---
310   gMC->Gstpar(idtmed[704],"LOSS",3.);
311   gMC->Gstpar(idtmed[704],"DRAY",1.);
312   // --- and in aluminium parts ---
313   gMC->Gstpar(idtmed[701],"LOSS",3.);
314   gMC->Gstpar(idtmed[701],"DRAY",1.);
315
316 }
317
318 //////////////////////////////////////////////////////////////////////////////
319
320 void AliPHOSv2::CreateGeometry()
321 {
322
323   AliPHOSv2 *PHOS_tmp = (AliPHOSv2*)gAlice->GetModule("PHOS");
324   if(PHOS_tmp==NULL){
325     
326     fprintf(stderr,"PHOS detector not found!\n");
327     return;
328     
329   }
330
331   // --- Dimensions of volumes ---
332   Float_t DPHOS[3], DPTXW[3], DPAIR[3];
333   Float_t DPUFP[3], DPUCP[3], DPASP[3], DPTIP[3], DPTXP[3];
334   Float_t DPTCB[3], DPCBL[3], DPSTC[3], DPPAP[3], DPXTL[3], DPSUP[3],
335     DPPIN[3];
336   Float_t DPCPV[3], DPCPA[3];
337
338   Float_t R, YO, XP1, YP1, PPHI, angle;
339   Int_t IDROTM[99];
340   Int_t i;
341   
342   Double_t const RADDEG=180.0/kPI;
343   //  Double_t const DEGRAD=kPI/180.0;
344
345   // --- Dimensions of PbWO4 crystal ---
346   //      PARAMETER(XTL_X=2.2,XTL_Y=18.,XTL_Z=2.2)
347   Float_t XTL_X=GetCrystalSize(0);
348   Float_t XTL_Y=GetCrystalSize(1);
349   Float_t XTL_Z=GetCrystalSize(2);
350
351   // --- Tyvek wrapper thickness
352   //      PARAMETER(PAP_THICK=0.01)
353   Float_t PAP_THICK=GetWrapThickness();
354
355   // --- Steel (titanium) cover thickness ---
356   //      PARAMETER(STE_THICK=0.005)
357   Float_t STE_THICK=GetPHOSextra(0);
358
359   // --- Crystal support height ---
360   //      PARAMETER(SUP_Y=6.95)
361   Float_t SUP_Y=GetPHOSextra(1);
362
363   // --- PIN-diode dimensions ---
364   //      PARAMETER(PIN_X=1.4,PIN_Y=0.4,PIN_Z=1.4)
365   Float_t PIN_X=GetPINSize(0);
366   Float_t PIN_Y=GetPINSize(1);
367   Float_t PIN_Z=GetPINSize(2);
368
369   // --- CPV thickness ---
370   //      PARAMETER(CPV_Y=0.5)
371   Float_t CPV_Y=GetCPVThickness();
372
373   // --- Foam Thermo Insulating outer cover dimensions ---
374   //      PARAMETER(FTI_X=214.6,FTI_Y=80.,FTI_Z=260.)
375   Float_t FTI_X=GetPHOSFoam(0);
376   Float_t FTI_Y=GetPHOSFoam(1);
377   Float_t FTI_Z=GetPHOSFoam(2);
378
379   // --- Thermo Insulating outer cover Upper plate thickness ---
380   //      PARAMETER(FTIU_THICK=4.)
381   Float_t FTIU_THICK=GetPHOSextra(2);
382
383   // --- Textolit Wall box dimentions ---
384   //      PARAMETER(TXW_X=209.,TXW_Y=71.,TXW_Z=250.)
385   Float_t TXW_X=GetPHOStxwall(0);
386   Float_t TXW_Y=GetPHOStxwall(1);
387   Float_t TXW_Z=GetPHOStxwall(2);
388
389   // --- Inner AIR volume dimensions ---
390   //      PARAMETER(AIR_X=206.,AIR_Y=66.,AIR_Z=244.)
391   Float_t AIR_X=GetPHOSAir(0);
392   Float_t AIR_Y=GetPHOSAir(1);
393   Float_t AIR_Z=GetPHOSAir(2);
394
395   // --- Upper Polystyrene Foam plate thickness ---
396   //      PARAMETER(UFP_Y=5.)
397   Float_t UFP_Y=GetPHOSextra(3);
398
399   // --- Thermo insulating Crystal Block wall thickness ---
400   //      PARAMETER(TCB_THICK=2.)
401   Float_t TCB_THICK=GetPHOSextra(4);
402
403   // --- Upper Cooling Plate thickness ---
404   //      PARAMETER(UCP_Y=0.06)
405   Float_t UCP_Y=GetPHOSextra(5);
406
407   // --- Al Support Plate thickness ---
408   //      PARAMETER(ASP_Y=10.)
409   Float_t ASP_Y=GetPHOSextra(6);
410
411   // --- Lower Thermo Insulating Plate thickness ---
412   //      PARAMETER(TIP_Y=3.)
413   Float_t TIP_Y=GetPHOSextra(7);
414
415   // --- Lower Textolit Plate thickness ---
416   //      PARAMETER(TXP_Y=1.)
417   Float_t TXP_Y=GetPHOSextra(8);
418
419   // --- 1/2 total gap between adjacent crystals
420   Float_t TOTAL_GAP=GetPHOSextra(9);
421
422   // --- Distance from IP to Foam Thermo Insulating top plate ---
423   //      PARAMETER(FTI_R=467.)
424   Float_t FTI_R=GetRadius(0);
425
426   // --- Distance from IP to Crystal Block top Surface (needs to be 460.) ---
427   //      PARAMETER(CBS_R=480.)
428   Float_t CBS_R=GetRadius(1);
429
430   // Get pointer to the array containing media indeces
431   Int_t *IDTMED = fIdtmed->GetArray()-699;
432
433   // --- Define PHOS box volume, fill with thermo insulating foam ---
434   DPHOS[0]=FTI_X/2.0;
435   DPHOS[1]=FTI_Y/2.0;
436   DPHOS[2]=FTI_Z/2.0;
437   gMC->Gsvolu("PHOS", "BOX ", IDTMED[706], DPHOS, 3);
438
439   // --- Define Textolit Wall box, position inside PHOS ---
440   DPTXW[0]=TXW_X/2.0;
441   DPTXW[1]=TXW_Y/2.0;
442   DPTXW[2]=TXW_Z/2.0;
443   gMC->Gsvolu("PTXW", "BOX ", IDTMED[707], DPTXW, 3);
444   YO=(FTI_Y-TXW_Y)/2.0-FTIU_THICK;
445   gMC->Gspos("PTXW", 1, "PHOS", 0.0, YO, 0.0, 0, "ONLY");
446
447   // --- Define Upper Polystyrene Foam Plate, place inside PTXW ---
448   // --- immediately below Foam Thermo Insulation Upper plate ---
449   DPUFP[0]=TXW_X/2.0;
450   DPUFP[1]=UFP_Y/2.0;
451   DPUFP[2]=TXW_Z/2.0;
452   gMC->Gsvolu("PUFP", "BOX ", IDTMED[703], DPUFP, 3);
453   YO=(TXW_Y-UFP_Y)/2.0;
454   gMC->Gspos("PUFP", 1, "PTXW", 0.0, YO, 0.0, 0, "ONLY");
455
456   // --- Define air-filled box, place inside PTXW ---
457   DPAIR[0]=AIR_X/2.0;
458   DPAIR[1]=AIR_Y/2.0;
459   DPAIR[2]=AIR_Z/2.0;
460   gMC->Gsvolu("PAIR", "BOX ", IDTMED[798], DPAIR, 3);
461   YO=(TXW_Y-AIR_Y)/2.0-UFP_Y;
462   gMC->Gspos("PAIR", 1, "PTXW", 0.0, YO, 0.0, 0, "ONLY");
463
464   // --- Define Thermo insulating Crystal Box, position inside PAIR ---
465   DPTCB[0]=GetNphi()*(XTL_X+2*TOTAL_GAP)/2.0+TCB_THICK;
466   DPTCB[1]=(XTL_Y+SUP_Y+PAP_THICK+STE_THICK)/2.0+TCB_THICK/2.0;
467   DPTCB[2]=GetNz()*(XTL_Z+2*TOTAL_GAP)/2.0+TCB_THICK;
468   gMC->Gsvolu("PTCB", "BOX ", IDTMED[706], DPTCB, 3);
469   YO=AIR_Y/2.0-DPTCB[1]-
470     (CBS_R-FTI_R-TCB_THICK-FTIU_THICK-UFP_Y);
471   gMC->Gspos("PTCB", 1, "PAIR", 0.0, YO, 0.0, 0, "ONLY");
472
473   // --- Define Crystal BLock filled with air, position it inside PTCB ---
474   DPCBL[0]=GetNphi()*(XTL_X+2*TOTAL_GAP)/2.0;
475   DPCBL[1]=(XTL_Y+SUP_Y+PAP_THICK+STE_THICK)/2.0;
476   DPCBL[2]=GetNz()*(XTL_Z+2*TOTAL_GAP)/2.0;
477   gMC->Gsvolu("PCBL", "BOX ", IDTMED[798], DPCBL, 3);
478   
479   // --- Divide PCBL in X (phi) and Z directions --
480   gMC->Gsdvn("PROW", "PCBL", GetNphi(), 1);
481   gMC->Gsdvn("PCEL", "PROW", GetNz(), 3);
482   YO=-TCB_THICK/2.0;
483   gMC->Gspos("PCBL", 1, "PTCB", 0.0, YO, 0.0, 0, "ONLY");
484
485   // --- Define STeel (actually, it's titanium) Cover volume, place inside PCEL
486   DPSTC[0]=(XTL_X+2*PAP_THICK)/2.0;
487   DPSTC[1]=(XTL_Y+SUP_Y+PAP_THICK+STE_THICK)/2.0;
488   DPSTC[2]=(XTL_Z+2*PAP_THICK+2*STE_THICK)/2.0;
489   gMC->Gsvolu("PSTC", "BOX ", IDTMED[704], DPSTC, 3);
490   gMC->Gspos("PSTC", 1, "PCEL", 0.0, 0.0, 0.0, 0, "ONLY");
491
492   // --- Define Tyvek volume, place inside PSTC ---
493   DPPAP[0]=XTL_X/2.0+PAP_THICK;
494   DPPAP[1]=(XTL_Y+SUP_Y+PAP_THICK)/2.0;
495   DPPAP[2]=XTL_Z/2.0+PAP_THICK;
496   gMC->Gsvolu("PPAP", "BOX ", IDTMED[702], DPPAP, 3);
497   YO=(XTL_Y+SUP_Y+PAP_THICK)/2.0-(XTL_Y+SUP_Y+PAP_THICK+STE_THICK)/2.0;
498   gMC->Gspos("PPAP", 1, "PSTC", 0.0, YO, 0.0, 0, "ONLY");
499
500   // --- Define PbWO4 crystal volume, place inside PPAP ---
501   DPXTL[0]=XTL_X/2.0;
502   DPXTL[1]=XTL_Y/2.0;
503   DPXTL[2]=XTL_Z/2.0;
504   gMC->Gsvolu("PXTL", "BOX ", IDTMED[699], DPXTL, 3);
505   YO=(XTL_Y+SUP_Y+PAP_THICK)/2.0-XTL_Y/2.0-PAP_THICK;
506   gMC->Gspos("PXTL", 1, "PPAP", 0.0, YO, 0.0, 0, "ONLY");
507
508   // --- Define crystal support volume, place inside PPAP ---
509   DPSUP[0]=XTL_X/2.0+PAP_THICK;
510   DPSUP[1]=SUP_Y/2.0;
511   DPSUP[2]=XTL_Z/2.0+PAP_THICK;
512   gMC->Gsvolu("PSUP", "BOX ", IDTMED[798], DPSUP, 3);
513   YO=SUP_Y/2.0-(XTL_Y+SUP_Y+PAP_THICK)/2.0;
514   gMC->Gspos("PSUP", 1, "PPAP", 0.0, YO, 0.0, 0, "ONLY");
515
516   // --- Define PIN-diode volume and position it inside crystal support ---
517   // --- right behind PbWO4 crystal
518   DPPIN[0]=PIN_X/2.0;
519   DPPIN[1]=PIN_Y/2.0;
520   DPPIN[2]=PIN_Z/2.0;
521   gMC->Gsvolu("PPIN", "BOX ", IDTMED[705], DPPIN, 3);
522   YO=SUP_Y/2.0-PIN_Y/2.0;
523   gMC->Gspos("PPIN", 1, "PSUP", 0.0, YO, 0.0, 0, "ONLY");
524
525   // --- Define Upper Cooling Panel, place it on top of PTCB ---
526   DPUCP[0]=DPTCB[0];
527   DPUCP[1]=UCP_Y/2.0;
528   DPUCP[2]=DPTCB[2];
529   gMC->Gsvolu("PUCP", "BOX ", IDTMED[701], DPUCP,3);
530   YO=(AIR_Y-UCP_Y)/2.0-(CBS_R-FTI_R-TCB_THICK-FTIU_THICK-UFP_Y-UCP_Y);
531   gMC->Gspos("PUCP", 1, "PAIR", 0.0, YO, 0.0, 0, "ONLY");
532
533   // --- Define Al Support Plate, position it inside PAIR ---
534   // --- right beneath PTCB ---
535   DPASP[0]=AIR_X/2.0;
536   DPASP[1]=ASP_Y/2.0;
537   DPASP[2]=AIR_Z/2.0;
538   gMC->Gsvolu("PASP", "BOX ", IDTMED[701], DPASP, 3);
539   YO=(AIR_Y-ASP_Y)/2.0-(CBS_R-FTI_R-FTIU_THICK-UFP_Y+DPCBL[1]*2);
540   gMC->Gspos("PASP", 1, "PAIR", 0.0, YO, 0.0, 0, "ONLY");
541
542   // --- Define Thermo Insulating Plate, position it inside PAIR ---
543   // --- right beneath PASP ---
544   DPTIP[0]=AIR_X/2.0;
545   DPTIP[1]=TIP_Y/2.0;
546   DPTIP[2]=AIR_Z/2.0;
547   gMC->Gsvolu("PTIP", "BOX ", IDTMED[706], DPTIP, 3);
548   YO=(AIR_Y-TIP_Y)/2.0-(CBS_R-FTI_R-FTIU_THICK-UFP_Y+DPCBL[1]*2+ASP_Y);
549   gMC->Gspos("PTIP", 1, "PAIR", 0.0, YO, 0.0, 0, "ONLY");
550
551   // --- Define Textolit Plate, position it inside PAIR ---
552   // --- right beneath PTIP ---
553   DPTXP[0]=AIR_X/2.0;
554   DPTXP[1]=TXP_Y/2.0;
555   DPTXP[2]=AIR_Z/2.0;
556   gMC->Gsvolu("PTXP", "BOX ", IDTMED[707], DPTXP, 3);
557   YO=(AIR_Y-TXP_Y)/2.0-
558     (CBS_R-FTI_R-FTIU_THICK-UFP_Y+DPCBL[1]*2+ASP_Y+TIP_Y);
559   gMC->Gspos("PTXP", 1, "PAIR", 0.0, YO, 0.0, 0, "ONLY");
560
561   // --- Define CPV volume, DON'T PLACE IT YET ---
562   // --- Divide in X and Z direction (same way as PCBL) ---
563   DPCPV[0]=DPCBL[0];
564   DPCPV[1]=CPV_Y/2.0;
565   DPCPV[2]=DPCBL[2];
566   //  gMC->Gsvolu("PCPV", "BOX ", IDTMED[700], DPCPV, 3);
567   gMC->Gsvolu("PCPV", "BOX ", IDTMED[798], DPCPV, 3);
568   gMC->Gsdvn("PCRO", "PCPV", GetNphi(), 1);
569   gMC->Gsdvn("PCCE", "PCRO", GetNz(), 3);
570
571   // Define CPV sensitive pad. It has the same size as PCCE.
572   DPCPA[0]=DPCBL[0]/GetNphi();
573   DPCPA[1]=CPV_Y/2.0;
574   DPCPA[2]=DPCBL[2]/GetNz();
575   gMC->Gsvolu("PCPA", "BOX ", IDTMED[700], DPCPA, 3);
576   gMC->Gspos("PCPA", 1, "PCCE", 0.0, 0.0, 0.0, 0, "ONLY");
577
578   // --- Position various PHOS units in ALICE setup ---
579   // --- PHOS itself first ---
580   PPHI=TMath::ATan(FTI_X/(2.0*FTI_R));
581   PPHI*=RADDEG;
582
583   for(i=1; i<=GetNModules(); i++){
584
585     angle=PPHI*2*(i-GetNModules()/2.0-0.5);
586     AliMatrix(IDROTM[i-1], 90.0, angle, 90.0, 90.0+angle, 0.0, 0.0);
587
588     // --- Position various PHOS units in ALICE setup ---
589     // --- PHOS itself first ---
590     R=FTI_R+FTI_Y/2.0;
591     XP1=R*TMath::Sin(angle/RADDEG);
592     YP1=-R*TMath::Cos(angle/RADDEG);
593     gMC->Gspos("PHOS", i, "ALIC", XP1, YP1, 0.0, IDROTM[i-1], "ONLY");
594
595     // --- Now position PCPV so that its plates are right on top of ---
596     // --- corresponding PHOS modules (previously called cradles) ---
597     R=FTI_R-CPV_Y/2.0;
598     XP1=R*TMath::Sin(angle/RADDEG);
599     YP1=-R*TMath::Cos(angle/RADDEG);
600     gMC->Gspos("PCPV", i, "ALIC", XP1, YP1, 0.0, IDROTM[i-1], "ONLY");
601     GetModuleAngle(i-1)=angle-90.0;
602
603   }
604
605   // --- Set volumes seen without their descendants for drawing ---
606   gMC->Gsatt("PCEL", "SEEN", -2);
607   gMC->Gsatt("PCCE", "SEEN", -2);
608
609 }
610
611 //////////////////////////////////////////////////////////////////////////////
612
613 void AliPHOSv2::StepManager(void)
614 {
615
616   Int_t blrc[4]; // (box, layer, row, column) indices
617   Float_t xyze[4]; // position wrt MRS and energy deposited
618   TLorentzVector pos;
619
620   Int_t *IDTMED=fIdtmed->GetArray()-699;
621
622   if(gMC->GetMedium()==IDTMED[700]){ // We are inside a CPV sensitive pad
623
624     gMC->TrackPosition(pos);
625     xyze[0]=pos[0];
626     xyze[1]=pos[1];
627     xyze[2]=pos[2];
628     xyze[3]=gMC->Edep();
629     
630     gMC->CurrentVolOffID(3, blrc[0]);
631     blrc[1]=1; // CPV corresponds to layer 1
632     gMC->CurrentVolOffID(2, blrc[2]);
633     gMC->CurrentVolOffID(1, blrc[3]);
634     
635     AddHit(gAlice->CurrentTrack(), blrc, xyze);
636
637   }
638
639   if(gMC->GetMedium()==IDTMED[699]){ // We are inside a PWO crystal
640
641     gMC->TrackPosition(pos);
642     xyze[0]=pos[0];
643     xyze[1]=pos[1];
644     xyze[2]=pos[2];
645     xyze[3]=gMC->Edep();
646
647     gMC->CurrentVolOffID(9, blrc[0]);
648     blrc[1]=2; // PWO crystals correspond to layer 2
649     gMC->CurrentVolOffID(4, blrc[2]);
650     gMC->CurrentVolOffID(3, blrc[3]);
651
652     AddHit(gAlice->CurrentTrack(), blrc, xyze);
653
654   }
655 }
656
657 ////////////////////////////////////////////////////////////////////////////
658
659 void AliPHOSv2::AddHit(Int_t track, Int_t *vol, Float_t *hits)
660 {
661
662   Int_t hitCounter;
663   TClonesArray &lhits = *fHits;
664   AliPHOShitv2 *newHit,*curHit;
665   
666   newHit=new AliPHOShitv2(fIshunt, track, vol, hits);
667   
668   for(hitCounter=0;hitCounter<fNhits;hitCounter++){
669     curHit=(AliPHOShitv2*)lhits[hitCounter];
670     if(*curHit==*newHit){
671       *curHit=*curHit+*newHit;
672       delete newHit;
673       return;
674     }
675   }
676   
677   new(lhits[fNhits++]) AliPHOShitv2(*newHit);
678   delete newHit;
679
680 }
681
682 ///////////////////////////////////////////////////////////////////////////////
683
684 ClassImp(AliPHOShitv2)
685
686 //////////////////////////////////////////////////////////////////////////////
687
688 AliPHOShitv2::AliPHOShitv2(Int_t shunt, Int_t track, Int_t *vol, Float_t *hits):
689 AliHit(shunt, track)
690 {
691
692    Int_t i;
693    for (i=0;i<4;i++)fVolume[i]=vol[i];
694
695    fX       = hits[0];
696    fY       = hits[1];
697    fZ       = hits[2];
698    fELOS    = hits[3];
699
700 }
701
702 //////////////////////////////////////////////////////////////////////////////
703
704 Bool_t AliPHOShitv2::operator==(AliPHOShitv2 const &rValue) const
705 {
706
707   Int_t volCounter;
708
709   //  if(fDet!=rValue.GetDet()) return kFALSE;
710   if(fTrack!=rValue.GetTrack()) return kFALSE;
711
712   for(volCounter=0;volCounter<4;volCounter++)
713     if(fVolume[volCounter]!=rValue.GetVolume(volCounter))return kFALSE;
714
715   return kTRUE;
716
717 }
718
719 /////////////////////////////////////////////////////////////////////////////
720
721 AliPHOShitv2 const AliPHOShitv2::operator+(AliPHOShitv2 const &rValue) const
722 {
723
724   AliPHOShitv2 added(*this);
725
726   added.fELOS+=rValue.GetEnergy();
727   return added;
728
729 }
730  
731 //////////////////////////////////////////////////////////////////////////////