]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHSegResV0.cxx
Avoid all warnings on SunOS. Again!
[u/mrichter/AliRoot.git] / RICH / AliRICHSegResV0.cxx
1 #include "AliRICHSegResV0.h"
2 #include "AliRun.h"
3 #include "TParticle.h"
4 #include "TMath.h"
5 #include "TRandom.h"
6
7
8 ClassImp(AliRICHsegmentation)
9 ClassImp(AliRICHresponse)       
10 //___________________________________________
11 ClassImp(AliRICHsegmentationV0)
12
13 void AliRICHsegmentationV0::Init(AliRICHchamber* Chamber)
14 {
15     fNpx=(Int_t) (Chamber->ROuter()/fDpx+1);
16     fNpy=(Int_t) (Chamber->ROuter()/fDpy+1);
17 }
18
19
20 Float_t AliRICHsegmentationV0::GetAnod(Float_t xhit)
21 {
22     Float_t wire= (xhit>0)? Int_t(xhit/fWireD)+0.5:Int_t(xhit/fWireD)-0.5;
23     return fWireD*wire;
24 }
25
26 void AliRICHsegmentationV0::SetPADSIZ(Float_t p1, Float_t p2)
27 {
28     fDpx=p1;
29     fDpy=p2;
30 }
31 void AliRICHsegmentationV0::GetPadIxy(Float_t x, Float_t y, Int_t &ix, Int_t &iy)
32 {
33 //  returns pad coordinates (ix,iy) for given real coordinates (x,y)
34 //
35     ix = (x>0)? Int_t(x/fDpx)+1 : Int_t(x/fDpx);
36     iy = (y>0)? Int_t(y/fDpy)+1 : Int_t(y/fDpy);
37     if (iy >  fNpy) iy= fNpy;
38     if (iy < -fNpy) iy=-fNpy;
39     if (ix >  fNpx) ix= fNpx;
40     if (ix < -fNpx) ix=-fNpx;
41 }
42 void AliRICHsegmentationV0::
43 GetPadCxy(Int_t ix, Int_t iy, Float_t &x, Float_t &y)
44 {
45 //  returns real coordinates (x,y) for given pad coordinates (ix,iy)
46 //
47     x = (ix>0) ? Float_t(ix*fDpx)-fDpx/2. : Float_t(ix*fDpx)-fDpx/2.;
48     y = (iy>0) ? Float_t(iy*fDpy)-fDpy/2. : Float_t(iy*fDpy)-fDpy/2.;
49 }
50
51 void AliRICHsegmentationV0::
52 FirstPad(Float_t xhit, Float_t yhit, Float_t dx, Float_t dy)
53 {
54     //
55     // Find the wire position (center of charge distribution)
56     Float_t x0a=GetAnod(xhit);
57     //
58     // and take fNsigma*sigma around this center
59     Float_t x01=x0a  - dx;
60     Float_t x02=x0a  + dx;
61     Float_t y01=yhit - dy;
62     Float_t y02=yhit + dy;
63     //
64     // find the pads over which the charge distributes
65     GetPadIxy(x01,y01,fixmin,fiymin);
66     GetPadIxy(x02,y02,fixmax,fiymax);    
67     // 
68     // Set current pad to lower left corner
69     fix=fixmin;
70     fiy=fiymin;
71     GetPadCxy(fix,fiy,fx,fy);
72 }
73
74 void AliRICHsegmentationV0::NextPad()
75 {
76     // 
77     // Step to next pad in integration region
78     if (fix != fixmax) {
79         fix++;
80     } else if (fiy != fiymax) {
81         fix=fixmin;
82         fiy++;
83     } else {
84         printf("\n Error: Stepping outside integration region\n ");
85     }
86     GetPadCxy(fix,fiy,fx,fy);
87 }
88
89 Int_t AliRICHsegmentationV0::MorePads()
90     
91 //
92 // Are there more pads in the integration region
93 {
94     if (fix == fixmax && fiy == fiymax) {
95         return 0;
96     } else {
97         return 1;
98         
99     }
100 }
101
102 void AliRICHsegmentationV0::SigGenInit(Float_t x,Float_t y,Float_t)
103 {
104 //
105 //  Initialises pad and wire position during stepping
106     fxt =x;
107     fyt =y;
108     GetPadIxy(x,y,fixt,fiyt);
109     fiwt=Int_t(x/fWireD)+1;
110     
111 }
112
113 Int_t AliRICHsegmentationV0::SigGenCond(Float_t x,Float_t y,Float_t)
114 {
115 //
116 //  Signal will be generated if particle crosses pad boundary or
117 //  boundary between two wires. 
118     Int_t ixt, iyt;
119     GetPadIxy(x,y,ixt,iyt);
120     Int_t iwt=Int_t(x/fWireD)+1;
121     
122     if ((ixt != fixt) || (iyt !=fiyt) || (iwt != fiwt)) {
123         return 1;
124     } else {
125         return 0;
126     }
127 }
128 void AliRICHsegmentationV0::
129 IntegrationLimits(Float_t& x1,Float_t& x2,Float_t& y1, Float_t& y2)
130 {
131     x1=fxt-fx-fDpx/2.;
132     x2=x1+fDpx;
133     y1=fyt-fy-fDpy/2.;
134     y2=y1+fDpy;    
135     
136 }
137
138 void AliRICHsegmentationV0::
139 Neighbours(Int_t iX, Int_t iY, Int_t* Nlist, Int_t Xlist[7], Int_t Ylist[7])
140 {
141 //Is used for the cluster finder, include diagonal elements
142     
143     *Nlist=4;Xlist[0]=Xlist[1]=iX;Xlist[2]=iX-1;Xlist[3]=iX+1;
144     Ylist[0]=iY-1;Ylist[1]=iY+1;Ylist[2]=Ylist[3]=iY;
145 }
146
147 void AliRICHsegmentationV0::
148 FitXY(AliRICHRecCluster* ,TClonesArray* )
149     // Default : Centre of gravity method
150 {
151     ;
152 }
153
154
155 //___________________________________________
156 ClassImp(AliRICHresponseV0)     
157     Float_t AliRICHresponseV0::IntPH(Float_t eloss)
158 {
159     // Get number of electrons and return charge
160     
161     Int_t nel;
162 //E9/26=magic number should parameter
163     nel= Int_t(eloss*1.e9/26.);
164     Float_t charge=0;
165     if (nel == 0) nel=1;
166     for (Int_t i=1;i<=nel;i++) {
167         charge -= fChslope*TMath::Log(gRandom->Rndm());    
168     }
169     return charge;
170 }
171 // -------------------------------------------
172 Float_t AliRICHresponseV0::IntXY(AliRICHsegmentation * segmentation)
173 {
174     
175     const Float_t invpitch = 1/fPitch;
176     Float_t response;
177 //
178 //  Integration limits defined by segmentation model
179 //  
180     
181     Float_t xi1, xi2, yi1, yi2;
182     segmentation->IntegrationLimits(xi1,xi2,yi1,yi2);
183     xi1=xi1*invpitch;
184     xi2=xi2*invpitch;
185     yi1=yi1*invpitch;
186     yi2=yi2*invpitch;
187     
188     //
189 // The Mathieson function 
190     Double_t ux1=fSqrtKx3*TMath::TanH(fKx2*xi1);
191     Double_t ux2=fSqrtKx3*TMath::TanH(fKx2*xi2);
192     
193     Double_t uy1=fSqrtKy3*TMath::TanH(fKy2*yi1);
194     Double_t uy2=fSqrtKy3*TMath::TanH(fKy2*yi2);
195     
196     response=4.*fKx4*(TMath::ATan(ux2)-TMath::ATan(ux1))*fKy4*(TMath::ATan(uy2)-TMath::ATan(uy1));
197
198     return response;       
199     
200 }
201 //___________________________________________
202 Int_t AliRICHresponseV0::FeedBackPhotons(Float_t source[3], Float_t qtot)
203 {
204   //
205   // Generate FeedBack photons
206   //
207   Int_t j, ipart, nt;
208     
209   //Probability of feedback
210   Float_t fAlphaFeed=0.05;
211     
212   Int_t sNfeed=0;
213   
214   // Local variables 
215   Float_t cthf, ranf[2], phif, enfp = 0, sthf, weight;
216   Int_t i, ifeed;
217   Float_t e1[3], e2[3], e3[3];
218   Float_t vmod, uswop;
219   Float_t fp, random;
220   Float_t dir[3], phi;
221   Int_t nfp;
222   Float_t pol[3], mom[3];
223   TLorentzVector position;
224   //
225   // Determine number of feedback photons
226
227   //  Get weight of current particle
228   TParticle *current = (TParticle*) 
229     (*gAlice->Particles())[gAlice->CurrentTrack()];
230     
231   ifeed = Int_t(current->GetWeight()/100+0.5);
232   ipart = gMC->TrackPid();
233   fp = fAlphaFeed * qtot;
234   nfp = gRandom->Poisson(fp);
235   
236   // This call to fill the time of flight
237   gMC->TrackPosition(position);
238   //
239   // Generate photons
240   for (i = 0; i <nfp; ++i) {
241         
242     // Direction
243     gMC->Rndm(ranf, 2);
244     cthf = ranf[0] * 2 - 1.;
245     if (cthf < 0)  continue;
246     sthf = TMath::Sqrt((1 - cthf) * (1 + cthf));
247     phif = ranf[1] * 2 * TMath::Pi();
248     //
249     gMC->Rndm(&random, 1);
250     if (random <= .57) {
251       enfp = 7.5e-9;
252     } else if (random <= .7) {
253       enfp = 6.4e-9;
254     } else {
255       enfp = 7.9e-9;
256     }
257
258     dir[0] = sthf * TMath::Sin(phif);
259     dir[1] = cthf;
260     dir[2] = sthf * TMath::Cos(phif);
261     gMC->Gdtom(dir, mom, 2);
262     mom[0]*=enfp;
263     mom[1]*=enfp;
264     mom[2]*=enfp;
265     
266     // Polarisation
267     e1[0] = 0;
268     e1[1] = -dir[2];
269     e1[2] = dir[1];
270     
271     e2[0] = -dir[1];
272     e2[1] = dir[0];
273     e2[2] = 0;
274     
275     e3[0] = dir[1];
276     e3[1] = 0;
277     e3[2] = -dir[0];
278     
279     vmod=0;
280     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
281     if (!vmod) for(j=0;j<3;j++) {
282       uswop=e1[j];
283       e1[j]=e3[j];
284       e3[j]=uswop;
285     }
286     vmod=0;
287     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
288     if (!vmod) for(j=0;j<3;j++) {
289       uswop=e2[j];
290       e2[j]=e3[j];
291       e3[j]=uswop;
292     }
293     
294     vmod=0;
295     for(j=0;j<3;j++) vmod+=e1[j]*e1[j];
296     vmod=TMath::Sqrt(1/vmod);
297     for(j=0;j<3;j++) e1[j]*=vmod;
298     
299     vmod=0;
300     for(j=0;j<3;j++) vmod+=e2[j]*e2[j];
301     vmod=TMath::Sqrt(1/vmod);
302     for(j=0;j<3;j++) e2[j]*=vmod;
303     
304     gMC->Rndm(ranf, 1);
305     phi = ranf[0] * 2 * TMath::Pi();
306     for(j=0;j<3;j++) pol[j]=e1[j]*TMath::Sin(phi)+e2[j]*TMath::Cos(phi);
307     gMC->Gdtom(pol, pol, 2);
308     
309     // Put photon on the stack and label it as feedback (51, 52) 
310     ++sNfeed;
311     if (ipart == 50000050 && ifeed != 50000052) {
312       weight = 5000;
313     } else {
314       weight = 5000;
315     }
316     gAlice->SetTrack(1, gAlice->CurrentTrack(), gMC->PDGFromId(50),
317                      mom,source,pol,position[3],
318                      "Cherenkov", nt, weight);
319   }
320   return(sNfeed);
321 }