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