]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv1.cxx
Minor corrections thanks to I.Hrivnacova
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv1.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 /////////////////////////////////////////////////////////
21 //  Manager and hits classes for set:PHOS version 1    //
22 /////////////////////////////////////////////////////////
23  
24 // --- ROOT system ---
25 #include "TH1.h"
26 #include "TRandom.h"
27 #include "TFile.h"
28 #include "TTree.h"
29 #include "TMath.h"
30
31 // --- galice header files ---
32 #include "AliPHOSv1.h"
33 #include "AliRun.h"
34
35 ClassImp(AliPHOSv1)
36
37 //______________________________________________________________________________
38
39
40 AliPHOSv1::AliPHOSv1()
41 {
42 }
43  
44 //______________________________________________________________________________
45
46 AliPHOSv1::AliPHOSv1(const char *name, const char *title)
47           : AliPHOS(name, title)
48 {
49 }
50  
51 //___________________________________________
52 void AliPHOSv1::CreateGeometry()
53 {
54
55   AliPHOS *PHOS_tmp = (AliPHOS*)gAlice->GetModule("PHOS");
56   if( NULL==PHOS_tmp )
57   {
58     printf("There isn't PHOS detector!\n");
59     return;
60   }
61 //  AliPHOS &PHOS = *PHOS_tmp;
62
63   //////////////////////////////////////////////////////////////////////////////
64
65   Int_t                 rotation_matrix_number=0;
66   Float_t               par[11],
67                         x,y,z;
68
69   const float           cell_length             = GetCrystalLength()+GetAirThickness()+GetWrapThickness()+GetPIN_Length(),
70                         cell_side_size          = GetCrystalSideSize()+2*GetAirThickness()+2*GetWrapThickness(),
71                         cradle_thikness         = cell_length;
72
73   //////////////////////////////////////////////////////////////////////////////
74   // CELL volume and subvolumes creation
75   //////////////////////////////////////////////////////////////////////////////
76
77   par[0] = GetCrystalSideSize()/2 + GetWrapThickness();
78   par[1] = GetCrystalSideSize()/2 + GetWrapThickness();
79   par[2] = GetCrystalLength()  /2 + GetWrapThickness()/2;
80   gMC->Gsvolu("WRAP","BOX ",GetPHOS_IDTMED_Tyvek(),par,3);
81
82   par[0] = GetCrystalSideSize()/2;
83   par[1] = GetCrystalSideSize()/2;
84   par[2] = GetCrystalLength()/2;
85   gMC->Gsvolu("CRST","BOX ",GetPHOS_IDTMED_PbWO4(),par,3);
86
87   // PIN
88   par[0] = GetPIN_SideSize()/2;
89   par[1] = GetPIN_SideSize()/2;
90   par[2] = GetPIN_Length()/2;
91   gMC->Gsvolu("PIN ","BOX ",GetPHOS_IDTMED_PIN(),par,3);
92
93   //////////////////////////////////////////////////////////////////////////////
94   // CRADLE creation.
95   //////////////////////////////////////////////////////////////////////////////
96
97   par[0] = cell_side_size/2 * GetNz();
98   par[1] = cell_side_size/2 * GetNphi();
99   par[2] = cradle_thikness/2;
100   gMC->Gsvolu("PHOS","BOX ",GetPHOS_IDTMED_AIR(),par,3);
101
102
103   par[0] = cell_side_size/2 * GetNz();
104   par[1] = cell_side_size/2 * GetNphi();
105   par[2] = cell_length/2;
106   gMC->Gsvolu("CRS0","BOX ",GetPHOS_IDTMED_AIR(),par,3);
107
108   x = 0;
109   y = 0;
110   z = -(cradle_thikness-cell_length)/2;
111   gMC->Gspos("CRS0",1,"PHOS",x,y,z,0,"ONLY");
112
113   gMC->Gsdvn("CRS1","CRS0",GetNphi(),2);
114   gMC->Gsdvn("CELL","CRS1",GetNz()  ,1);
115
116   //////////////////////////////////////////////////////////////////////////////
117   // CELL creation
118   //////////////////////////////////////////////////////////////////////////////
119
120   x = 0;
121   y = 0;
122   z = -GetWrapThickness()/2;
123   gMC->Gspos("CRST",1,"WRAP",x,y,z,0,"ONLY");
124
125   x = 0;
126   y = 0;
127   z = GetPIN_Length()/2;
128   gMC->Gspos("WRAP",1,"CELL",x,y,z,0,"ONLY");
129
130   x = 0;
131   y = 0;
132   z = -GetCrystalLength()/2-GetWrapThickness()/2;
133   gMC->Gspos("PIN ",1,"CELL",x,y,z,0,"ONLY");
134
135   //////////////////////////////////////////////////////////////////////////////
136   // CELL has been created.
137   //////////////////////////////////////////////////////////////////////////////
138
139
140   //////////////////////////////////////////////////////////////////////////////
141   // End of CRADLE creation.
142   //////////////////////////////////////////////////////////////////////////////
143
144
145   //////////////////////////////////////////////////////////////////////////////
146   // PHOS creation
147   //////////////////////////////////////////////////////////////////////////////
148
149   for( int i=0; i<GetCradlesAmount(); i++ )
150   {
151     Float_t cradle_angle     = 27.,
152             cradle_angle_pos = -90+(i-(GetCradlesAmount()-1)/2.) *
153                                (cradle_angle+GetAngleBetweenCradles());
154     // Cradles are numerated in clock reversed order. (general way of angle increment)
155
156     Float_t r = GetRadius() + cradle_thikness/2;
157     x = r*cos(cradle_angle_pos*kPI/180);
158     y = r*sin(cradle_angle_pos*kPI/180);
159     z = 0;
160     AliMatrix(rotation_matrix_number, 0,0 , 90,90+cradle_angle_pos , 90,180+cradle_angle_pos);
161     gMC->Gspos("PHOS",i+1,"ALIC",x,y,z,rotation_matrix_number,"ONLY");
162
163     GetCradleAngle(i) = cradle_angle_pos;
164   }
165   AddPHOSCradles();
166
167   //////////////////////////////////////////////////////////////////////////////
168   // All is done.
169   // Print some information.
170   //////////////////////////////////////////////////////////////////////////////
171 }
172
173 void AliPHOSv1::StepManager()
174 {
175   static Bool_t inwold=0;   // Status of previous ctrak->inwvol
176   Int_t copy;
177
178   int cradle_number, cell_Z, cell_Phi;  // Variables that describe cell position.
179
180   if( gMC->GetMedium() == GetPHOS_IDTMED_PIN() && (gMC->IsTrackInside() || gMC->IsTrackExiting()==2) && inwold && gMC->TrackCharge()!=0 )
181   {
182     // GEANT particle just have entered into PIN diode.
183
184     AliPHOS &PHOS = *(AliPHOS*)gAlice->GetModule("PHOS");
185
186     gMC->CurrentVolOffID(4,copy);
187     cradle_number  = copy-1;
188     gMC->CurrentVolOffID(1,copy);
189     cell_Z         = copy-1;
190     gMC->CurrentVolOffID(2,copy);
191     cell_Phi       = copy-1;
192 /*
193         cradle_number  = cvolu->number[cvolu->nlevel-5]-1;
194         cell_Z         = cvolu->number[cvolu->nlevel-2]-1;
195         cell_Phi       = cvolu->number[cvolu->nlevel-3]-1;
196 */
197
198     TH2S &h = PHOS.GetCradle(cradle_number).fChargedTracksInPIN;
199     h.AddBinContent(h.GetBin(cell_Z,cell_Phi));
200   }
201
202   //////////////////////////////////////////////////////////////////////////////
203
204   if( gMC->GetMedium() == GetPHOS_IDTMED_PbWO4() )
205   {
206     // GEANT particle into crystal.
207
208     AliPHOS &PHOS = *(AliPHOS*)gAlice->GetModule("PHOS");
209
210     gMC->CurrentVolOffID(5,copy);
211     cradle_number  = copy-1;
212     gMC->CurrentVolOffID(2,copy);
213     cell_Z         = copy-1;
214     gMC->CurrentVolOffID(3,copy);
215     cell_Phi       = copy-1;
216 /*
217         cradle_number  = cvolu->number[cvolu->nlevel-6]-1;
218         cell_Z         = cvolu->number[cvolu->nlevel-3]-1;
219         cell_Phi       = cvolu->number[cvolu->nlevel-4]-1;
220 */
221     TH2F &h = PHOS.GetCradle(cradle_number).fCellEnergy;
222     h.AddBinContent(h.GetBin(cell_Z,cell_Phi),gMC->Edep());
223   }
224
225   //////////////////////////////////////////////////////////////////////////////
226
227
228   inwold=gMC->IsTrackEntering();         // Save current status of GEANT variable.
229 }
230