]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSv1.cxx
Eliminate useless include
[u/mrichter/AliRoot.git] / PHOS / AliPHOSv1.cxx
1 /////////////////////////////////////////////////////////
2 //  Manager and hits classes for set:PHOS version 1    //
3 /////////////////////////////////////////////////////////
4  
5 // --- ROOT system ---
6 #include "TH1.h"
7 #include "TRandom.h"
8 #include "TFile.h"
9 #include "TTree.h"
10 #include "TMath.h"
11
12 // --- galice header files ---
13 #include "AliPHOSv1.h"
14 #include "AliRun.h"
15
16 ClassImp(AliPHOSv1)
17
18 //______________________________________________________________________________
19
20
21 AliPHOSv1::AliPHOSv1()
22 {
23 }
24  
25 //______________________________________________________________________________
26
27 AliPHOSv1::AliPHOSv1(const char *name, const char *title)
28           : AliPHOS(name, title)
29 {
30 }
31  
32 //___________________________________________
33 void AliPHOSv1::CreateGeometry()
34 {
35
36   AliPHOS *PHOS_tmp = (AliPHOS*)gAlice->GetModule("PHOS");
37   if( NULL==PHOS_tmp )
38   {
39     printf("There isn't PHOS detector!\n");
40     return;
41   }
42 //  AliPHOS &PHOS = *PHOS_tmp;
43
44   //////////////////////////////////////////////////////////////////////////////
45
46   Int_t                 rotation_matrix_number=0;
47   Float_t               par[11],
48                         x,y,z;
49
50   const float           cell_length             = GetCrystalLength()+GetAirThickness()+GetWrapThickness()+GetPIN_Length(),
51                         cell_side_size          = GetCrystalSideSize()+2*GetAirThickness()+2*GetWrapThickness(),
52 //                        cell_angle              = 180/kPI * 2 * atan(cell_side_size/2 / GetRadius()),        // radians
53                         cradle_thikness         = cell_length + GetCPV_Thickness() + GetCPV_PHOS_Distance(),
54                         distance_to_CPV         = GetRadius() - GetCPV_Thickness() - GetCPV_PHOS_Distance();
55
56   //////////////////////////////////////////////////////////////////////////////
57   // CELL volume and subvolumes creation
58   //////////////////////////////////////////////////////////////////////////////
59
60   par[0] = GetCrystalSideSize()/2 + GetWrapThickness();
61   par[1] = GetCrystalSideSize()/2 + GetWrapThickness();
62   par[2] = GetCrystalLength()  /2 + GetWrapThickness()/2;
63   gMC->Gsvolu("WRAP","BOX ",GetPHOS_IDTMED_Tyvek(),par,3);
64
65   par[0] = GetCrystalSideSize()/2;
66   par[1] = GetCrystalSideSize()/2;
67   par[2] = GetCrystalLength()/2;
68   gMC->Gsvolu("CRST","BOX ",GetPHOS_IDTMED_PbWO4(),par,3);
69
70   // PIN
71   par[0] = GetPIN_SideSize()/2;
72   par[1] = GetPIN_SideSize()/2;
73   par[2] = GetPIN_Length()/2;
74   gMC->Gsvolu("PIN ","BOX ",GetPHOS_IDTMED_PIN(),par,3);
75
76   //////////////////////////////////////////////////////////////////////////////
77   // CRADLE,CPV creation.
78   //////////////////////////////////////////////////////////////////////////////
79
80   par[0] = cell_side_size/2 * GetNz();
81   par[1] = cell_side_size/2 * GetNphi();
82   par[2] = cradle_thikness/2;
83   gMC->Gsvolu("PHOS","BOX ",GetPHOS_IDTMED_AIR(),par,3);
84
85 //par[0] : the same as above
86 //par[1] : the same as above
87   par[2] = GetCPV_Thickness()/2;
88   gMC->Gsvolu("CPV ","BOX ",GetPHOS_IDTMED_CPV(),par,3);
89
90   x = 0;
91   y = 0;
92   z = (cell_length+GetCPV_PHOS_Distance())/2;
93   gMC->Gspos("CPV ",1,"PHOS",x,y,z,0,"ONLY");
94
95   par[0] = cell_side_size/2 * GetNz();
96   par[1] = cell_side_size/2 * GetNphi();
97   par[2] = cell_length/2;
98   gMC->Gsvolu("CRS0","BOX ",GetPHOS_IDTMED_AIR(),par,3);
99
100   x = 0;
101   y = 0;
102   z = -(cradle_thikness-cell_length)/2;
103   gMC->Gspos("CRS0",1,"PHOS",x,y,z,0,"ONLY");
104
105   gMC->Gsdvn("CRS1","CRS0",GetNphi(),2);
106   gMC->Gsdvn("CELL","CRS1",GetNz()  ,1);
107
108   //////////////////////////////////////////////////////////////////////////////
109   // CELL creation
110   //////////////////////////////////////////////////////////////////////////////
111
112   x = 0;
113   y = 0;
114   z = -GetWrapThickness()/2;
115   gMC->Gspos("CRST",1,"WRAP",x,y,z,0,"ONLY");
116
117   x = 0;
118   y = 0;
119   z = GetPIN_Length()/2;
120   gMC->Gspos("WRAP",1,"CELL",x,y,z,0,"ONLY");
121
122   x = 0;
123   y = 0;
124   z = -GetCrystalLength()/2-GetWrapThickness()/2;
125   gMC->Gspos("PIN ",1,"CELL",x,y,z,0,"ONLY");
126
127   //////////////////////////////////////////////////////////////////////////////
128   // CELL has been created.
129   //////////////////////////////////////////////////////////////////////////////
130
131 //   int n=0;
132 //   z = -(GetCPV_Thickness()+GetCPV_PHOS_Distance())/2;
133 //
134 //   for( int iy=0; iy<GetNphi(); iy++ )
135 //   {
136 //     y = (iy-(GetNphi()-1)/2.)*cell_side_size;
137 //     for( int ix=0; ix<GetNz(); ix++ )
138 //     {
139 //       x = (ix-(GetNz()-1)/2.)*cell_side_size;
140 //       gMC->Gspos("CELL",++n,"PHOS",x,y,z,0,"ONLY");
141 //     }
142 //   }
143
144   //////////////////////////////////////////////////////////////////////////////
145   // End of CRADLE creation.
146   //////////////////////////////////////////////////////////////////////////////
147
148
149   //////////////////////////////////////////////////////////////////////////////
150   // PHOS creation
151   //////////////////////////////////////////////////////////////////////////////
152
153   for( int i=0; i<GetCradlesAmount(); i++ )
154   {
155     float c                = distance_to_CPV,           // Distance to CPV
156           l                = cell_side_size*GetNphi()/2,      // Cradle half size around beam (for rect. geom.)
157           cradle_angle     = 360/kPI*atan(l/c),
158           cradle_angle_pos = -90+(i-(GetCradlesAmount()-1)/2.) * (cradle_angle+GetAngleBetweenCradles());
159     // Cradles are numerated in clock reversed order. (general way of angle increment)
160
161     float   r       = GetRadius() + cradle_thikness/2;
162     x = r*cos(cradle_angle_pos*kPI/180);
163     y = r*sin(cradle_angle_pos*kPI/180);
164     z = 0;
165     AliMatrix(rotation_matrix_number, 0,0 , 90,90+cradle_angle_pos , 90,180+cradle_angle_pos);
166     gMC->Gspos("PHOS",i+1,"ALIC",x,y,z,rotation_matrix_number,"ONLY");
167
168     GetCradleAngle(i) = cradle_angle_pos;
169 //
170 //    int n = PHOS.fCradles->GetEntries();
171 //    PHOS.fCradles->Add(new AliPHOSCradle( 1,            // geometry.
172 //                                          GetCrystalSideSize    (),
173 //                                          GetCrystalLength      (),
174 //                                          GetWrapThickness      (),
175 //                                          GetAirThickness       (),
176 //                                          GetPIN_SideSize       (),
177 //                                          GetPIN_Length         (),
178 //                                          GetRadius             (),
179 //                                          GetCPV_Thickness      (),
180 //                                          GetCPV_PHOS_Distance  (),
181 //                                          GetNz                 (),
182 //                                          GetNphi               (),
183 //                                          cradle_angle_pos      ));
184 //
185 //    if( n+1 != PHOS.fCradles->GetEntries() ||
186 //        NULL == PHOS.fCradles->At(n) )
187 //    {
188 //      cout << "  Can not create or add AliPHOSCradle.\n";
189 //      exit(1);
190 //    }
191   }
192   AddPHOSCradles();
193
194   //////////////////////////////////////////////////////////////////////////////
195   // All is done.
196   // Print some information.
197   //////////////////////////////////////////////////////////////////////////////
198 }
199
200 void AliPHOSv1::StepManager()
201 {
202   static Bool_t inwold=0;   // Status of previous ctrak->inwvol
203   Int_t copy;
204
205   int cradle_number, cell_Z, cell_Phi;  // Variables that describe cell position.
206
207   if( gMC->GetMedium() == GetPHOS_IDTMED_PIN() && (gMC->IsTrackInside() || gMC->IsTrackExiting()==2) && inwold && gMC->TrackCharge()!=0 )
208   {
209     // GEANT particle just have entered into PIN diode.
210
211     AliPHOS &PHOS = *(AliPHOS*)gAlice->GetModule("PHOS");
212
213     gMC->CurrentVolOffID(4,copy);
214     cradle_number  = copy-1;
215     gMC->CurrentVolOffID(1,copy);
216     cell_Z         = copy-1;
217     gMC->CurrentVolOffID(2,copy);
218     cell_Phi       = copy-1;
219 /*
220         cradle_number  = cvolu->number[cvolu->nlevel-5]-1;
221         cell_Z         = cvolu->number[cvolu->nlevel-2]-1;
222         cell_Phi       = cvolu->number[cvolu->nlevel-3]-1;
223 */
224
225     TH2S &h = PHOS.GetCradle(cradle_number).fChargedTracksInPIN;
226     h.AddBinContent(h.GetBin(cell_Z,cell_Phi));
227   }
228
229   //////////////////////////////////////////////////////////////////////////////
230
231   if( gMC->GetMedium() == GetPHOS_IDTMED_PbWO4() )
232   {
233     // GEANT particle into crystal.
234
235     AliPHOS &PHOS = *(AliPHOS*)gAlice->GetModule("PHOS");
236
237     gMC->CurrentVolOffID(5,copy);
238     cradle_number  = copy-1;
239     gMC->CurrentVolOffID(2,copy);
240     cell_Z         = copy-1;
241     gMC->CurrentVolOffID(3,copy);
242     cell_Phi       = copy-1;
243 /*
244         cradle_number  = cvolu->number[cvolu->nlevel-6]-1;
245         cell_Z         = cvolu->number[cvolu->nlevel-3]-1;
246         cell_Phi       = cvolu->number[cvolu->nlevel-4]-1;
247 */
248     TH2F &h = PHOS.GetCradle(cradle_number).fCellEnergy;
249     h.AddBinContent(h.GetBin(cell_Z,cell_Phi),gMC->Edep());
250   }
251
252   //////////////////////////////////////////////////////////////////////////////
253
254   if( gMC->GetMedium()==GetPHOS_IDTMED_CPV() && (gMC->IsTrackInside() || gMC->IsTrackExiting()) && inwold )
255   {
256     // GEANT particle just have entered into CPV detector.
257
258     AliPHOS &PHOS = *(AliPHOS*)gAlice->GetModule("PHOS");
259
260     gMC->CurrentVolOffID(1,cradle_number);
261     cradle_number--;
262 //        cradle_number  = cvolu->number[cvolu->nlevel-2]-1;
263
264     // Save CPV x,y hits position of charged particles.
265
266     AliPHOSCradle  &cradle = PHOS.GetCradle(cradle_number);
267
268     TLorentzVector xyz;
269     TVector3 v;
270     gMC->TrackPosition(xyz);
271
272     float x,y,l;
273     float R = cradle.GetRadius() - cradle.GetCPV_PHOS_Distance() - cradle.GetCPV_Thikness();
274     cradle.GetXY(xyz.Vect(),v,R,x,y,l);
275     if( PHOS.fDebugLevel>0 )
276       if( l<0 )
277         printf("PHOS_STEP:  warning: negative distance to CPV!! %f\n", l);
278
279     // Store current particle in the list of Cradle particles.
280     TLorentzVector  pmom;
281     gMC->TrackMomentum(pmom);
282     float     Px      =       pmom[0],
283               Py      =       pmom[1],
284               Pz      =       pmom[2];
285     Int_t     Ipart   =       gMC->TrackPid();
286
287 //     TClonesArray &P=cradle.GetParticles();
288 //     new( P[P.GetEntries()] ) AliPHOSgamma(x,0,y,0,ctrak->getot,0,Px,Py,Pz);
289     cradle.GetParticles().Add(new AliPHOSgamma(x,y,gMC->Etot(),Px,Py,Pz,Ipart));
290
291     if( gMC->TrackCharge()!=0 )
292       cradle.AddCPVHit(x,y);
293   }
294
295   inwold=gMC->IsTrackEntering();         // Save current status of GEANT variable.
296 }
297