]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/Opticals.h
fgMCEvGen changed to fMCEvGen;
[u/mrichter/AliRoot.git] / RICH / Opticals.h
1 #ifdef __CINT__
2 void Opticals()
3 {
4   gROOT->Reset();
5 #endif   
6   int i;
7   const Int_t kNbins=26;
8   Float_t aPckov[kNbins];
9   for(i=0;i<kNbins;i++){    //Photons energy intervals
10     aPckov[i]=(Float_t(i)*0.1+5.5)*1e-9;
11   }
12   
13   Float_t aIndexFreon[kNbins];  
14   Float_t aIndexQuartz[kNbins];
15   Float_t aIndexOpaqueQuartz[kNbins];
16   Float_t aIndexCH4[kNbins];
17   Float_t aIndexGrid[kNbins];
18         
19   Float_t  e1= 10.666;Float_t  e2= 18.125;  Float_t  f1= 46.411; Float_t  f2= 228.71;
20   for (i=0;i<kNbins;i++){
21     aIndexFreon[i]        = aPckov[i] * .0172 * 1e9 + 1.177;
22     Float_t ene=aPckov[i]*1e9;
23     Float_t a=f1/(e1*e1 - ene*ene);
24     Float_t b=f2/(e2*e2 - ene*ene);
25     aIndexQuartz[i]        = TMath::Sqrt(1. + a + b );
26     aIndexOpaqueQuartz[i]  =1;
27     aIndexCH4[i]      =1.000444;
28     aIndexGrid[i]         =1;
29   } 
30     
31   Float_t aAbsFreon[kNbins]={179.0987,  179.0987,    179.0987,  179.0987,   179.0987,  
32                              179.0987,  179.0987,    179.0987,  179.0987,   142.9206, 
33                               56.64957,  25.58622,    13.95293,  12.03905,   10.42953, 
34                                8.804196,  7.069031,    4.461292,  2.028366,   1.293013, 
35                                0.577267,  0.40746,     0.334964,  0.0,        0.0,       
36                                0};
37     
38
39   Float_t aAbsQuartz[kNbins]={105.8,   65.52,  48.58,  42.85,   35.79,
40                                31.262, 28.598, 27.527, 25.007,  22.815, 
41                                21.004, 19.266, 17.525, 15.878,  14.177, 
42                                11.719, 9.282,   6.62,   4.0925,  2.601, 
43                                1.149,  0.667,   0.3627, 0.192,   0.1497, 
44                                0.10857};
45     
46   Float_t aAbsOpaqueQuartz[kNbins];
47   Float_t aAbsCH4[kNbins];
48   Float_t aAbsGrid[kNbins];
49   Float_t aAbsCsI[kNbins];
50     for (i=0;i<kNbins;i++){
51       aAbsCH4[i]    =AbsoCH4(aPckov[i]*1e9); 
52       aAbsOpaqueQuartz[i]=1e-5; 
53       aAbsCsI[i]        =1e-4; 
54       aAbsGrid[i]       =1e-4; 
55     }
56     
57   Float_t aQeCsI[kNbins] = {0.000199999995, 0.000600000028, 0.000699999975, 0.00499999989, 0.00749999983,
58                             0.010125,       0.0242999997,   0.0405000001,   0.0688500032,  0.105299994, 
59                             0.121500008,    0.141749993,    0.157949999,    0.162,         0.166050002, 
60                             0.167669997,    0.174299985,    0.176789999,    0.179279998,   0.182599992, 
61                             0.18592,        0.187579989,    0.189239994,    0.190899998,   0.207499996, 
62                             0.215799987};
63   Float_t aQeAll[kNbins];
64   for(i=0;i<kNbins;i++){
65     aQeCsI[i]/= (1.0-Fresnel(aPckov[i]*1e9,1.0,0)); //FRESNEL LOSS CORRECTION
66     aQeAll[i]=1; //QE for all other materials except for PC must be 1.
67   }
68        
69                     
70 #ifdef __CINT__
71
72   
73 //Freon, Quartz, Opaque Quartz,Methane,CsI,Grid   
74   const Int_t kFreonMarker= 24; const Int_t kFreonColor=kRed;  
75   const Int_t kCH4Marker= 25;   const Int_t kCH4Color=kGreen;  
76   const Int_t kSiO2Marker= 26;  const Int_t kSiO2Color=kBlue;  
77   const Int_t kCsIMarker = 2;   const Int_t kCsIColor =kMagenta;  
78   
79   TCanvas *pC=new TCanvas("c1","RICH optics to check",800,900);
80
81   pC->Divide(2,2);           
82            
83   pC->cd(1);                 
84   TGraph *pAbsFreonGr=new TGraph(kNbins,aPckov,aAbsFreon);
85   pAbsFreonGr->SetMarkerStyle(kFreonMarker); pAbsFreonGr->SetMarkerColor(kFreonColor);
86   TGraph *pAbsSiO2Gr=new TGraph(kNbins,aPckov,aAbsQuartz);
87   pAbsSiO2Gr->SetMarkerStyle(kSiO2Marker); pAbsSiO2Gr->SetMarkerColor(kSiO2Color);  
88   TMultiGraph *pAbsMG=new TMultiGraph();
89   TLegend *pAbsLegend=new TLegend(0.6,0.3,0.85,0.5);
90   pAbsMG->Add(pAbsFreonGr);     pAbsLegend->AddEntry(pAbsFreonGr,   "freon","p");            //1
91   pAbsMG->Add(pAbsSiO2Gr);      pAbsLegend->AddEntry(pAbsSiO2Gr,  "quartz","p");          //2
92   pAbsMG->Draw("APL");                                  
93   pAbsMG->GetXaxis()->SetTitle("energy, GeV");
94   pAbsMG->GetYaxis()->SetTitle("absorption length, cm");
95   pAbsMG->Draw("APL");
96   pAbsLegend->Draw();   
97
98
99   pC->cd(2);
100   TGraph *pIndexFreonGr=new TGraph(kNbins,aPckov,aIndexFreon);
101   pIndexFreonGr->SetMarkerStyle(kFreonMarker); pIndexFreonGr->SetMarkerColor(kFreonColor);  
102   TGraph *pIndexSiO2Gr=new TGraph(kNbins,aPckov,aIndexQuartz);
103   pIndexSiO2Gr->SetMarkerStyle(kSiO2Marker); pIndexSiO2Gr->SetMarkerColor(kSiO2Color);  
104   TGraph *pIndexCH4Gr=new TGraph(kNbins,aPckov,aIndexCH4);
105   pIndexCH4Gr->SetMarkerStyle(kCH4Marker); pIndexCH4Gr->SetMarkerColor(kCH4Color);   
106   TMultiGraph *pIndexMG=new TMultiGraph();
107   TLegend *pIndexLegend=new TLegend(0.6,0.2,0.85,0.4);
108   pIndexMG->Add(pIndexFreonGr);     pIndexLegend->AddEntry(pIndexFreonGr,   "freon","p");            //1
109   pIndexMG->Add(pIndexSiO2Gr);      pIndexLegend->AddEntry(pIndexSiO2Gr,  "quartz","p");          
110   pIndexMG->Add(pIndexCH4Gr);       pIndexLegend->AddEntry(pIndexCH4Gr,   "CH4","p");          
111   pIndexMG->Draw("APL");                                
112   pIndexMG->GetXaxis()->SetTitle("energy, GeV");
113   pIndexMG->GetYaxis()->SetTitle("refraction index");
114   pIndexMG->Draw("APL");
115   pIndexLegend->Draw();      
116
117   pC->cd(3);      
118   TGraph *pAbsCH4Gr=new TGraph(kNbins,aPckov,aAbsCH4);
119   pAbsCH4Gr->SetMarkerStyle(kCH4Marker); pAbsCH4Gr->SetMarkerColor(kCH4Color);
120   pAbsCH4Gr->Draw("APL");
121   pAbsCH4Gr->GetXaxis()->SetTitle("energy, GeV");
122   pAbsCH4Gr->SetTitle("CH4 absorption length, cm");
123   pAbsCH4Gr->Draw("APL");
124
125   
126   pC->cd(4);  
127   TGraph *pQeCsIGr=new TGraph(kNbins,aPckov,aQeCsI);
128   pQeCsIGr->SetMarkerStyle(kCsIMarker); pQeCsIGr->SetMarkerColor(kCsIColor);  
129   pQeCsIGr->Draw("APL");
130   pQeCsIGr->GetXaxis()->SetTitle("energy, GeV");
131   pQeCsIGr->SetTitle("CsI QE, all others are 1");
132   pQeCsIGr->Draw("APL");
133 }//main
134
135 Float_t AbsoCH4(Float_t x)
136 {
137
138     //KLOSCH,SCH4(9),WL(9),EM(9),ALENGTH(31)
139     Float_t sch4[9] = {.12,.16,.23,.38,.86,2.8,7.9,28.,80.};              //MB X 10^22
140     //Float_t wl[9] = {153.,152.,151.,150.,149.,148.,147.,146.,145};
141     Float_t em[9] = {8.1,8.158,8.212,8.267,8.322,8.378,8.435,8.493,8.55};
142     const Float_t kLosch=2.686763E19;                                      // LOSCHMIDT NUMBER IN CM-3
143     const Float_t kIgas1=100, kIgas2=0, kOxy=10., kWater=5., kPressure=750.,kTemperature=283.;                                      
144     Float_t pn=kPressure/760.;
145     Float_t tn=kTemperature/273.16;
146     
147         
148 // ------- METHANE CROSS SECTION -----------------
149 // ASTROPH. J. 214, L47 (1978)
150         
151     Float_t sm=0;
152     if (x<7.75) 
153         sm=.06e-22;
154     
155     if(x>=7.75 && x<=8.1)
156     {
157         Float_t c0=-1.655279e-1;
158         Float_t c1=6.307392e-2;
159         Float_t c2=-8.011441e-3;
160         Float_t c3=3.392126e-4;
161         sm=(c0+c1*x+c2*x*x+c3*x*x*x)*1.e-18;
162     }
163     
164     if (x> 8.1)
165     {
166         Int_t j=0;
167         while (x<=em[j] && x>=em[j+1])
168         {
169             j++;
170             Float_t a=(sch4[j+1]-sch4[j])/(em[j+1]-em[j]);
171             sm=(sch4[j]+a*(x-em[j]))*1e-22;
172         }
173     }
174     
175     Float_t dm=(kIgas1/100.)*(1.-((kOxy+kWater)/1.e6))*kLosch*pn/tn;
176     Float_t abslm=1./sm/dm;
177     
178 //    ------- ISOBUTHANE CROSS SECTION --------------
179 //     i-C4H10 (ai) abs. length from curves in
180 //     Lu-McDonald paper for BARI RICH workshop .
181 //     -----------------------------------------------------------
182     
183     Float_t ai;
184     Float_t absli;
185     if (kIgas2 != 0) 
186     {
187         if (x<7.25)
188             ai=100000000.;
189         
190         if(x>=7.25 && x<7.375)
191             ai=24.3;
192         
193         if(x>=7.375)
194             ai=.0000000001;
195         
196         Float_t si = 1./(ai*kLosch*273.16/293.);                    // ISOB. CRO.SEC.IN CM2
197         Float_t di=(kIgas2/100.)*(1.-((kOxy+kWater)/1.e6))*kLosch*pn/tn;
198         absli =1./si/di;
199     }
200     else
201         absli=1.e18;
202 //    ---------------------------------------------------------
203 //
204 //       transmission of O2
205 //
206 //       y= path in cm, x=energy in eV
207 //       so= cross section for UV absorption in cm2
208 //       do= O2 molecular density in cm-3
209 //    ---------------------------------------------------------
210     
211     Float_t abslo;
212     Float_t so=0;
213     if(x>=6.0)
214     {
215         if(x>=6.0 && x<6.5)
216         {
217             so=3.392709e-13 * TMath::Exp(2.864104 *x);
218             so=so*1e-18;
219         }
220         
221         if(x>=6.5 && x<7.0) 
222         {
223             so=2.910039e-34 * TMath::Exp(10.3337*x);
224             so=so*1e-18;
225         }
226             
227
228         if (x>=7.0) 
229         {
230             Float_t a0=-73770.76;
231             Float_t a1=46190.69;
232             Float_t a2=-11475.44;
233             Float_t a3=1412.611;
234             Float_t a4=-86.07027;
235             Float_t a5=2.074234;
236             so= a0+(a1*x)+(a2*x*x)+(a3*x*x*x)+(a4*x*x*x*x)+(a5*x*x*x*x*x);
237             so=so*1e-18;
238         }
239         
240         Float_t dox=(kOxy/1e6)*kLosch*pn/tn;
241         abslo=1./so/dox;
242     }
243     else
244         abslo=1.e18;
245 //     ---------------------------------------------------------
246 //
247 //       transmission of H2O
248 //
249 //       y= path in cm, x=energy in eV
250 //       sw= cross section for UV absorption in cm2
251 //       dw= H2O molecular density in cm-3
252 //     ---------------------------------------------------------
253     
254     Float_t abslw;
255     
256     Float_t b0=29231.65;
257     Float_t b1=-15807.74;
258     Float_t b2=3192.926;
259     Float_t b3=-285.4809;
260     Float_t b4=9.533944;
261     
262     if(x>6.75)
263     {    
264         Float_t sw= b0+(b1*x)+(b2*x*x)+(b3*x*x*x)+(b4*x*x*x*x);
265         sw=sw*1e-18;
266         Float_t dw=(kWater/1e6)*kLosch*pn/tn;
267         abslw=1./sw/dw;
268     }
269     else
270         abslw=1.e18;
271             
272 //    ---------------------------------------------------------
273     
274     Float_t alength=1./(1./abslm+1./absli+1./abslo+1./abslw);
275     return (alength);
276 }//AbsoCH4
277
278
279 Float_t Fresnel(Float_t ene,Float_t pdoti, Bool_t pola)
280 {
281
282     //ENE(EV), PDOTI=COS(INC.ANG.), PDOTR=COS(POL.PLANE ROT.ANG.)
283     
284     Float_t en[36] = {5.0,5.1,5.2,5.3,5.4,5.5,5.6,5.7,5.8,5.9,6.0,6.1,6.2,
285                       6.3,6.4,6.5,6.6,6.7,6.8,6.9,7.0,7.1,7.2,7.3,7.4,7.5,7.6,7.7,
286                       7.8,7.9,8.0,8.1,8.2,8.3,8.4,8.5};
287      
288
289     Float_t csin[36] = {2.14,2.21,2.33,2.48,2.76,2.97,2.99,2.59,2.81,3.05,
290                         2.86,2.53,2.55,2.66,2.79,2.96,3.18,3.05,2.84,2.81,2.38,2.11,
291                         2.01,2.13,2.39,2.73,3.08,3.15,2.95,2.73,2.56,2.41,2.12,1.95,
292                         1.72,1.53};
293       
294     Float_t csik[36] = {0.,0.,0.,0.,0.,0.196,0.408,0.208,0.118,0.49,0.784,0.543,
295                         0.424,0.404,0.371,0.514,0.922,1.102,1.139,1.376,1.461,1.253,0.878,
296                         0.69,0.612,0.649,0.824,1.347,1.571,1.678,1.763,1.857,1.824,1.824,
297                         1.714,1.498};
298     Float_t xe=ene;
299     Int_t  j=Int_t(xe*10)-49;
300     Float_t cn=csin[j]+((csin[j+1]-csin[j])/0.1)*(xe-en[j]);
301     Float_t ck=csik[j]+((csik[j+1]-csik[j])/0.1)*(xe-en[j]);
302
303     //FORMULAE FROM HANDBOOK OF OPTICS, 33.23 OR
304     //W.R. HUNTER, J.O.S.A. 54 (1964),15 , J.O.S.A. 55(1965),1197
305
306     Float_t sinin=TMath::Sqrt(1-pdoti*pdoti);
307     Float_t tanin=sinin/pdoti;
308
309     Float_t c1=cn*cn-ck*ck-sinin*sinin;
310     Float_t c2=4*cn*cn*ck*ck;
311     Float_t aO=TMath::Sqrt(0.5*(TMath::Sqrt(c1*c1+c2)+c1));
312     Float_t b2=0.5*(TMath::Sqrt(c1*c1+c2)-c1);
313     
314     Float_t rs=((aO-pdoti)*(aO-pdoti)+b2)/((aO+pdoti)*(aO+pdoti)+b2);
315     Float_t rp=rs*((aO-sinin*tanin)*(aO-sinin*tanin)+b2)/((aO+sinin*tanin)*(aO+sinin*tanin)+b2);
316     
317
318     //CORRECTION FACTOR FOR SURFACE ROUGHNESS
319     //B.J. STAGG  APPLIED OPTICS, 30(1991),4113
320
321     Float_t sigraf=18.;
322     Float_t lamb=1240/ene;
323     Float_t fresn;
324  
325     Float_t  rO=TMath::Exp(-(4*TMath::Pi()*pdoti*sigraf/lamb)*(4*TMath::Pi()*pdoti*sigraf/lamb));
326
327     if(pola)
328     {
329         Float_t pdotr=0.8;                                 //DEGREE OF POLARIZATION : 1->P , -1->S
330         fresn=0.5*(rp*(1+pdotr)+rs*(1-pdotr));
331     }
332     else
333         fresn=0.5*(rp+rs);
334       
335     fresn = fresn*rO;
336     return(fresn);
337 }//Fresnel(...)
338 #endif