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