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