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