]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHParam.cxx
Removing static data members (TF1 objects) from AliRICHParam to make it working on...
[u/mrichter/AliRoot.git] / RICH / AliRICHParam.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 #include "AliRICHParam.h" //class header
16 #include "AliESD.h"
17 #include <TCanvas.h>      //TestXXX() 
18 #include <TLatex.h>
19 #include <THStack.h>
20 #include <TLegend.h>
21 #include <TView.h>
22 #include <TPolyMarker3D.h>
23 #include <TPolyLine3D.h>
24 #include <TPolyLine.h>
25 #include <TSystem.h>
26 #include <TVector2.h>
27 #include <TVector3.h>
28 #include <TRotation.h>
29 #include <AliCDBManager.h> //CdbRead()
30 #include <AliCDBStorage.h> //CdbRead()
31 #include <AliCDBEntry.h>   //CdbRead()
32 #include <AliRunLoader.h>  //Stack()
33 #include <AliStack.h>      //Stack()
34 #include <TParticle.h>     //Stack()    
35 #include "AliRICHHelix.h"  //TestTrans()
36
37 ClassImp(AliRICHParam)
38 AliRICHParam * AliRICHParam::fgInstance             =0x0;     //singleton pointer               
39
40 Double_t AliRICHParam::fgMass[5]              ={0.00051,0.10566,0.13957,0.49360,0.93828};  
41
42
43 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
44 AliRICHParam::AliRICHParam():TNamed("RichParam","default version") 
45 {
46 // Here all the intitializition is taken place when AliRICHParam::Instance() is invoked for the first time.
47 // In particulare, matrices to be used for LORS<->MARS trasnformations are initialized from TGeo structure.    
48 // Note that TGeoManager should be already initialized from geometry.root file  
49   for(Int_t iCh=0;iCh<kNchambers;iCh++) fMatrix[iCh]=(TGeoHMatrix*)gGeoManager->GetVolume("ALIC")->GetNode(Form("RICH_%i",iCh+1))->GetMatrix();
50   CdbRead(0,0);
51   fgInstance=this; 
52 }//ctor
53 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
54 Float_t AliRICHParam::AbsCH4(Float_t eV)
55 {
56 // Evaluate the absorbtion lenght of CH4 for a photon of energy eV in electron-volts
57   const Float_t kLoschmidt=2.686763e19;                                      // LOSCHMIDT NUMBER IN CM-3
58   const Float_t kPressure=750.0;          //mm of Hg
59   const Float_t kTemperature=283.0;       //K (10 grad C)                               
60   const Float_t kPn=kPressure/760.;
61   const Float_t kTn=kTemperature/273.16;
62   const Float_t kC0=-1.655279e-1;
63   const Float_t kC1= 6.307392e-2;
64   const Float_t kC2=-8.011441e-3;
65   const Float_t kC3= 3.392126e-4;
66                 
67   Float_t crossSection=0;                        
68   if (eV<7.75) 
69     crossSection=0.06e-22;
70   else                 //------ METHANE CROSS SECTION cm-2 ASTROPH. J. 214, L47 (1978)                                               
71     crossSection=(kC0+kC1*eV+kC2*eV*eV+kC3*eV*eV*eV)*1.e-18;
72     
73     Float_t density=kLoschmidt*kPn/kTn; //CH4 molecular concentration (cm^-3)
74     return 1.0/(density*crossSection);
75 }//AbsoCH4()
76 //__________________________________________________________________________________________________sss
77 void AliRICHParam::CdbRead(Int_t run,Int_t version)
78 {
79 // This methode read all the calibration information and initialise corresponding fields for requested run number
80 // Arguments: run - run number for which to retrieve calibration
81 //            version- version number   
82 //   Returns: none      
83
84   AliCDBEntry *pEntry=AliCDBManager::Instance()->Get("RICH/RICHConfig/RefIdxC6F14",run,0,version); //try to get from common local storage  
85   if(pEntry){
86     fIdxC6F14=(TF2*)pEntry->GetObject(); delete pEntry;
87   }else{
88     AliWarning("No valid calibarion, the hardcoded will be used!");
89     fIdxC6F14=new TF2("RidxC4F14","sqrt(1+0.554*(1239.84e-9/x)^2/((1239.84e-9/x)^2-5796)-0.0005*(y-20))",5.5e-9,8.5e-9,0,50); //DiMauro mail
90     fIdxC6F14->SetUniqueID(20);//T=20 deg C
91   }
92 }//CdbRead()
93 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
94 void AliRICHParam::Print(Option_t* opt) const
95 {
96 // print some usefull (hopefully) info on some internal guts of RICH parametrisation 
97   Printf("Pads in chamber (%3i,%3i) in sector (%2i,%2i) pad size (%4.2f,%4.2f)",NpadsX(),NpadsY(),NpadsXsec(),NpadsYsec(),PadSizeX(),PadSizeY());
98   
99   for(Int_t i=0;i<kNchambers;i++) fMatrix[i]->Print(opt);
100 }//Print()
101 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
102 void AliRICHParam::TestSeg()
103 {
104 // Provides a set of pictures to test segementation currently in use.    
105 // Arguments: none
106 //   Returns: none    
107   new TCanvas("pads","PC segmentation - pads display",700,600);
108   gPad->Range(-5,-5,PcSizeX()+5,PcSizeY()+15);
109   TVector p(2);   TVector2 c;    TVector2 b;   //current: pad, pad center, pad boundary
110 // list of corners:
111   Double_t x0=0,x1=SecSizeX(),x2=SecSizeX()+DeadZone()                                                      ,x3=PcSizeX();  
112   Double_t y0=0,y1=SecSizeY(),y2=SecSizeY()+DeadZone(),y3=2*SecSizeY()+DeadZone(),y4=PcSizeY()-SecSizeY(),y5=PcSizeY();  
113   DrawSectors();
114 //header 
115   TLatex t;
116   t.SetTextSize(0.02); t.SetTextColor(kBlack); t.SetTextAlign(11);
117   t.DrawLatex(0,PcSizeY()+10,Form("IP in front of this page. pad size %.2fx%.2fcm   dead zone %.2fcm",PadSizeX(),PadSizeY(),DeadZone()));
118   t.DrawLatex(0,PcSizeY()+ 5,Form("Pc  %.2fx%.2f cm %ix%i pads               Sec %.2fx%.2f cm %ix%i pads",
119                                          PcSizeX()     , PcSizeY()     , NpadsX()    , NpadsY()                                 ,
120                                          SecSizeX() , SecSizeY() , NpadsXsec() , NpadsYsec()                              ));
121 //sectors  
122   t.SetTextSize(0.015); t.SetTextColor(kRed); t.SetTextAlign(22);
123   c=Pad2Loc( 40, 24); t.DrawText(c.X(),c.Y(),Form("sec 1 (%.2f,%.2f)",c.X(),c.Y()  ));  
124   c=Pad2Loc( 40, 75); t.DrawText(c.X(),c.Y(),Form("sec 3 (%.2f,%.2f)",c.X(),c.Y()  ));  
125   c=Pad2Loc( 40,121); t.DrawText(c.X(),c.Y(),Form("sec 5 (%.2f,%.2f)",c.X(),c.Y()  ));  
126   c=Pad2Loc(120, 24); t.DrawText(c.X(),c.Y(),Form("sec 2 (%.2f,%.2f)",c.X(),c.Y()  ));  
127   c=Pad2Loc(120, 75); t.DrawText(c.X(),c.Y(),Form("sec 4 (%.2f,%.2f)",c.X(),c.Y()  ));  
128   c=Pad2Loc(120,121); t.DrawText(c.X(),c.Y(),Form("sec 6 (%.2f,%.2f)",c.X(),c.Y()  ));  
129 //coners  
130   t.SetTextSize(0.015); t.SetTextColor(kBlue);
131
132   b.Set(x0,y0);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(11);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
133   b.Set(x0,y1);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(13);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
134   b.Set(x0,y2);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(11);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
135   b.Set(x0,y3);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(13);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
136   b.Set(x0,y4);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(11);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
137   b.Set(x0,y5);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(13);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
138   
139   b.Set(x1,y0);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(31);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
140   b.Set(x1,y1);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(33);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
141   b.Set(x1,y2);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(31);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
142   b.Set(x1,y3);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(33);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
143   b.Set(x1,y4);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(31);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
144   b.Set(x1,y5);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(33);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
145   
146   b.Set(x2,y0);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(11);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
147   b.Set(x2,y1);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(13);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
148   b.Set(x2,y2);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(11);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
149   b.Set(x2,y3);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(13);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
150   b.Set(x2,y4);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(11);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
151   b.Set(x2,y5);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(13);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
152   
153   b.Set(x3,y0);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(31);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
154   b.Set(x3,y1);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(33);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
155   b.Set(x3,y2);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(31);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
156   b.Set(x3,y3);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(33);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
157   b.Set(x3,y4);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(31);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
158   b.Set(x3,y5);p=Loc2Pad(b);c=Pad2Loc(p);t.SetTextAlign(33);t.DrawText(c.X(),c.Y(),Form("(%.2f,%.2f)-(%.0f,%.0f)-(%.2f,%.2f)",b.X(),b.Y(),p(0),p(1),c.X(),c.Y()));
159 }//TestSeg()
160 //__________________________________________________________________________________________________
161 void AliRICHParam::TestResp()
162 {
163 // Provides a set of plot to check the response parametrisation currently in use.  
164 // Arguments: none
165 //   Returns: none    
166   TCanvas *pC=new TCanvas("c","Amplification test",900,800);
167   pC->Divide(1,2);
168   
169   
170   const Int_t kNpoints=8;
171   THStack *pStackPhot=new THStack("StackPhot","photons");
172   THStack *pStackMip =new THStack("StackMip","mips");
173   TLegend *pLeg=new TLegend(0.6,0.2,0.9,0.5,"legend");    
174   TH1F *apHphot[kNpoints];
175   TH1F *apHmip[kNpoints];
176   
177   Double_t starty=0;
178   Double_t deltay=AliRICHParam::SecSizeY()/kNpoints;
179   
180   for(int i=0;i<kNpoints;i++){
181     apHphot[i]=new TH1F(Form("hphot%i",i),"Qdc for Photon;QDC;Counts",500,0,500); apHphot[i]->SetLineColor(i);pStackPhot->Add(apHphot[i]);                 
182     apHmip[i] =new TH1F(Form("hmip%i",i),"Qdc for Mip;QDC;Counts",4000,0,4000);   apHmip[i]->SetLineColor(i);pStackMip->Add(apHmip[i]);                 
183     
184     pLeg->AddEntry(apHphot[i],Form("@(10,%5.2f->%5.2f)",starty+i*deltay,starty+i*deltay-SecSizeY()/2));
185   }
186         
187   
188   TVector2 x2(0,0);  
189   for(Int_t i=0;i<10000;i++){//events loop
190     for(int j=0;j<kNpoints;j++){
191       x2.Set(10,starty+j*deltay);
192       apHphot[j]->Fill(TotQdc(x2,0));
193       apHmip[j]->Fill(TotQdc(x2,gRandom->Landau(600,150)*1e-9));
194     }
195   }
196   
197   pC->cd(1);  pStackMip->Draw("nostack");
198   pC->cd(2);  pStackPhot->Draw("nostack"); pLeg->Draw();
199 }//TestResp()
200 //__________________________________________________________________________________________________
201 void AliRICHParam::TestTrans()
202 {
203 // Tests transformation methods
204 // Arguments: none
205 //   Returns: none    
206   
207   AliRICHParam *pParam=AliRICHParam::Instance();
208   Int_t iNpointsX=50,iNpointsY=50;  
209   new TCanvas("trasform","Test LORS-MARS transform"); TLatex t; t.SetTextSize(0.02);
210   
211   TView *pView=new TView(1);  pView->SetRange(-400,-400,-400,400,400,400);
212   DrawAxis();  
213   for(Int_t iCham=1;iCham<=7;iCham++){//chamber loop
214     AliRICHHelix helix(2.5,Norm(iCham).Theta()*TMath::RadToDeg(),Norm(iCham).Phi()*TMath::RadToDeg());
215     helix.RichIntersect(AliRICHParam::Instance());
216     TPolyMarker3D *pChamber=new TPolyMarker3D(iNpointsX*iNpointsY);
217     Int_t i=0;
218     for(Double_t x=0;x<PcSizeX();x+=PcSizeX()/iNpointsX)
219       for(Double_t y=0;y<PcSizeY();y+=PcSizeY()/iNpointsY){//step loop
220         TVector3 v3=pParam->Lors2Mars(iCham,x,y,kPc); TVector2 v2=pParam->Mars2Lors(iCham,v3,kPc);//LORS->MARS->LORS
221         Double_t dx=v2.X()-x , dy=v2.Y()-y; 
222         if(dx>0.000001 || dy>0.000001) Printf("Problem in MARS<->LORS transformations dx=%f dy=%f!!!",dx,dy);
223         pChamber->SetPoint(i++,v3.X(),v3.Y(),v3.Z());//Pc plane point in MARS
224       }//step loop
225     pChamber->SetMarkerSize(1);
226     pChamber->SetMarkerColor(iCham);
227     pChamber->Draw();  
228     helix.Draw();
229     t.SetNDC();t.SetTextColor(iCham); t.DrawText(0.1,iCham*0.1,Form("Chamber %i",iCham));    
230   }//chambers loop   
231 }//TestTrans()
232 //__________________________________________________________________________________________________
233 void AliRICHParam::DrawAxis()
234 {
235 // This utility methode draws axis on geometry scene  
236 // Arguments: none
237 //   Returns: none    
238   Double_t x[6]={0,0,0,300,0,0};  Double_t y[6]={0,0,0,0,300,0};  Double_t z[6]={0,0,0,0,0,300};  
239   TPolyLine3D *pXaxis=new TPolyLine3D(2,x);pXaxis->SetLineColor(kRed);   pXaxis->Draw();
240   TPolyLine3D *pYaxis=new TPolyLine3D(2,y);pYaxis->SetLineColor(kGreen); pYaxis->Draw();
241   TPolyLine3D *pZaxis=new TPolyLine3D(2,z);pZaxis->SetLineColor(kBlue);  pZaxis->Draw();  
242 }
243 //__________________________________________________________________________________________________
244 void AliRICHParam::DrawSectors() 
245
246 // Utility methode draws RICH chamber sectors on event display.
247 // Arguments: none
248 //   Returns: none      
249   Double_t xLeft[5]  = {0,0,SecSizeX(),SecSizeX(),0};
250   Double_t xRight[5] = {SecSizeX()+DeadZone(),SecSizeX()+DeadZone(),PcSizeX(),PcSizeX(),SecSizeX()+DeadZone()};
251   
252   Double_t yDown[5]   = {0,SecSizeY(),SecSizeY(),0,0};
253   Double_t yCenter[5] = {  SecSizeY()+DeadZone(),2*SecSizeY()+DeadZone(),2*SecSizeY()+DeadZone(),
254                            SecSizeY()+DeadZone(),SecSizeY()+DeadZone()};  
255   Double_t yUp[5]     = {2*SecSizeY()+2*DeadZone(),PcSizeY(),PcSizeY(),2*SecSizeY()+2*DeadZone(),2*SecSizeY()+2*DeadZone()};
256   
257   TPolyLine *sec1 = new TPolyLine(5,xLeft ,yDown);    sec1->SetLineColor(21);  sec1->Draw();
258   TPolyLine *sec2 = new TPolyLine(5,xRight,yDown);    sec2->SetLineColor(21);  sec2->Draw();
259   TPolyLine *sec3 = new TPolyLine(5,xLeft ,yCenter);  sec3->SetLineColor(21);  sec3->Draw();
260   TPolyLine *sec4 = new TPolyLine(5,xRight,yCenter);  sec4->SetLineColor(21);  sec4->Draw();
261   TPolyLine *sec5 = new TPolyLine(5,xLeft, yUp);      sec5->SetLineColor(21);  sec5->Draw();
262   TPolyLine *sec6 = new TPolyLine(5,xRight,yUp);      sec6->SetLineColor(21);  sec6->Draw();
263 }//DrawSectors()
264 //__________________________________________________________________________________________________
265 TVector3 AliRICHParam::ForwardTracing(TVector3 entranceTrackPoint, TVector3 vectorTrack, Double_t thetaC, Double_t phiC)
266 {
267 // Trace a single Ckov photon from a given emission point up to photocathode taking into account ref indexes of materials it travereses
268   TVector3 vBad(-999,-999,-999);
269   TVector3 nPlane(0,0,1);
270   Double_t planeZposition = 0.5*RadThick();
271   TVector3 planePoint(0,0,0.5*RadThick()); //this is plane parallel to window which contains emission point 
272   TVector3 emissionPoint = PlaneIntersect(vectorTrack,entranceTrackPoint,nPlane,planePoint);
273   Double_t thetaout,phiout;
274   AnglesInDRS(vectorTrack.Theta(),vectorTrack.Phi(),thetaC,phiC,thetaout,phiout);
275   TVector3 vectorPhotonInC6F14;  
276   vectorPhotonInC6F14.SetMagThetaPhi(1,thetaout,phiout);
277   planeZposition=RadThick();
278   planePoint.SetXYZ(0,0,planeZposition);
279   TVector3 entranceToSiO2Point = PlaneIntersect(vectorPhotonInC6F14,emissionPoint,nPlane,planePoint);
280
281   Double_t photonEn = EckovMean();
282   Double_t angleInSiO2 = SnellAngle(IdxC6F14(EckovMean()),IdxSiO2(EckovMean()),vectorPhotonInC6F14.Theta());if(angleInSiO2<0) return vBad;
283   TVector3 vectorPhotonInSiO2;
284   vectorPhotonInSiO2.SetMagThetaPhi(1,angleInSiO2,phiout);
285 //  planeZposition+=AliRICHParam::SiO2Thickness();
286   planeZposition+=WinThick();
287   planePoint.SetXYZ(0,0,planeZposition);
288   TVector3 entranceToCH4 = PlaneIntersect(vectorPhotonInSiO2,entranceToSiO2Point,nPlane,planePoint);
289 //  entranceToCH4.Dump();
290
291   //  Double_t angleInCH4 = SnellAngle(AliRICHParam::IndOfRefSiO2(6.755),AliRICHParam::IndOfRefCH4,angleInSiO2);
292   Double_t angleInCH4 = SnellAngle(IdxSiO2(photonEn),IdxCH4(photonEn),vectorPhotonInSiO2.Theta());if(angleInCH4<0) return vBad;
293   TVector3 vectorPhotonInCH4;
294   vectorPhotonInCH4.SetMagThetaPhi(1,angleInCH4,phiout);
295 //  planeZposition+=AliRICHParam::GapProx();
296   planeZposition+=Pc2Win();
297   planePoint.SetXYZ(0,0,planeZposition);
298   TVector3 impactToPC = PlaneIntersect(vectorPhotonInCH4,entranceToCH4,nPlane,planePoint);
299 //  impactToPC.Dump();
300   return impactToPC;
301 }//FowardTracing
302 //__________________________________________________________________________________________________
303 TVector3 AliRICHParam::PlaneIntersect(const TVector3 &lineDir,const TVector3 &linePoint,const TVector3 &planeNorm,const TVector3 &planePoint)
304 {
305 // Finds an intersection point between a line and plane.
306 // Arguments:  lineDir,linePoint    - vector along the line and any point of the line
307 //             planeNorm,planePoint - vector normal to the plane and any point of the plane
308 //   Returns:  point of intersection if any
309   if(planeNorm*lineDir==0) return   TVector3(-999,-999,-999);
310   TVector3 diff=planePoint-linePoint;
311   Double_t sint=(planeNorm*diff)/(planeNorm*lineDir);
312   return linePoint+sint*lineDir;
313 }//PlaneIntersect
314 //__________________________________________________________________________________________________ 
315 Double_t AliRICHParam::SnellAngle(Float_t n1, Float_t n2, Float_t theta1)
316 {
317 // Compute the angle of refraction out of Snell law
318 // Arguments: n1 - ref idx of first substance  
319 //            n2 - ref idx of second substance  
320 //            n1 - photon impact angle in the first substance i.e. angle between the photon direction and vector normal to the surface (radians)
321 //   Returns: photon refraction angle, i.e. angle in the second substance (radians)
322   Double_t sinref=(n1/n2)*TMath::Sin(theta1);
323   if(sinref>1.)    return -999;  
324   else             return TMath::ASin(sinref);
325 }//SnellAngle
326 //__________________________________________________________________________________________________
327 void AliRICHParam::AnglesInDRS(Double_t trackTheta,Double_t trackPhi,Double_t thetaCerenkov,Double_t phiCerenkov,Double_t &tout,Double_t &pout)
328 {
329 // Setup the rotation matrix of the track...
330
331   TRotation mtheta;
332   TRotation mphi;
333   TRotation minv;
334   TRotation mrot;
335   
336   mtheta.RotateY(trackTheta);
337   mphi.RotateZ(trackPhi);
338   
339   mrot = mphi * mtheta;
340     //  minv = mrot.Inverse();
341
342   TVector3 photonInRadiator(1,1,1);
343
344   photonInRadiator.SetTheta(thetaCerenkov);
345   photonInRadiator.SetPhi(phiCerenkov);
346   photonInRadiator = mrot * photonInRadiator;
347   tout=photonInRadiator.Theta();
348   pout=photonInRadiator.Phi();
349 }//AnglesInDRS
350 /*
351 void DrawRing()
352 {
353
354   //  Float_t xGraph[1000],yGraph[1000];
355
356   Float_t type;
357   Float_t MassOfParticle;
358   Float_t beta;
359   Float_t nfreon;
360
361   Float_t ThetaCerenkov;
362
363   Float_t Xtoentr = GetEntranceX();
364   Float_t Ytoentr = GetEntranceY();
365
366   Float_t pmod = GetTrackMomentum();
367   Float_t TrackTheta = GetTrackTheta();
368   Float_t TrackPhi = GetTrackPhi();
369
370   SetPhotonEnergy(AliRICHParam::MeanCkovEnergy());
371   SetFreonRefractiveIndex();
372
373   SetEmissionPoint(RadiatorWidth/2.);
374
375   ThetaCerenkov = GetThetaCerenkov();
376   FindBetaFromTheta(ThetaCerenkov);
377   nfreon = GetFreonRefractiveIndex();
378   
379   Int_t nPoints = 100;
380
381   Int_t nPointsToDraw = 0;
382   for(Int_t i=0;i<nPoints;i++)
383     {
384       Float_t phpad = 2*TMath::Pi()*i/nPoints;
385       SetThetaPhotonInTRS(thetacer);
386       SetPhiPhotonInTRS(phpad);
387       FindPhotonAnglesInDRS();
388       Float_t Radius = FromEmissionToCathode();
389       if (Radius == 999.) continue;
390       xGraph[nPointsToDraw] = GetXPointOnCathode() + GetShiftX();
391       yGraph[nPointsToDraw] = GetYPointOnCathode() + GetShiftY();
392       nPointsToDraw++;
393     }
394   gra = new TGraph(nPointsToDraw,xGraph,yGraph);
395   gra->Draw("AC"); 
396 }
397 //__________________________________________________________________________________________________
398 */
399 void AliRICHParam::TestHit2SDigs(Double_t x,Double_t y,Double_t e,Bool_t isNew)
400 {
401 //Test  hit->sdigits procedures
402 //Arguments: isNew - if true use new (abs pad) procedure else use old one (TVector)
403 //  Returns: none
404   TClonesArray *pSDigLst=new TClonesArray("AliRICHDigit");
405   Int_t iQtot=-1;
406   if(isNew){
407     iQtot=Hit2SDigs(10101,e,pSDigLst);        //new technique
408   }else{
409     iQtot=Hit2SDigs(TVector2(x,y),e,pSDigLst);//old technique
410   }
411   pSDigLst->Print();
412   Double_t dQsum=0;
413   for(Int_t i=0;i<pSDigLst->GetEntriesFast();i++)
414     dQsum+=((AliRICHDigit*)pSDigLst->At(i))->Qdc();
415   Printf("Qtot=%i Qsum=%.2f ",iQtot,dQsum);
416 }
417 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
418 Int_t AliRICHParam::Stack(Int_t evt,Int_t tid)
419 {
420 // Prints some usefull info from stack
421 // Arguments: evt - event number. if not -1 print info only for that event
422 //            tid - track id. if not -1 then print it and all it's mothers if any   
423 //   Returns: mother tid of the given tid if any
424   AliRunLoader *pAL=AliRunLoader::Open(); 
425   if(pAL->LoadHeader()) return -1;
426   if(pAL->LoadKinematics()) return -1;
427   
428   Int_t mtid=-1;
429   Int_t iNevt=pAL->GetNumberOfEvents();     Printf("This session contains %i event(s)",iNevt);
430   
431   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
432     if(evt!=-1 && evt!=iEvt) continue; //in case one needs to print the requested event, ignore all others
433     pAL->GetEvent(iEvt);    
434     AliStack *pStack=pAL->Stack();  
435     if(tid==-1){                        //print all tids for this event
436       for(Int_t i=0;i<pStack->GetNtrack();i++) pStack->Particle(i)->Print();
437       Printf("totally %i tracks including %i primaries for event %i out of %i event(s)",pStack->GetNtrack(),pStack->GetNprimary(),iEvt,iNevt);
438     }else{                              //print only this tid and it;s mothers
439       if(tid<0 || tid>pStack->GetNtrack()) {Printf("Wrong tid, valid tid range for event %i is 0-%i",iEvt,pStack->GetNtrack());break;}
440       TParticle *pTrack=pStack->Particle(tid); mtid=pTrack->GetFirstMother();
441       TString str=pTrack->GetName();
442       while((tid=pTrack->GetFirstMother()) >= 0){
443         pTrack=pStack->Particle(tid);
444         str+=" from ";str+=pTrack->GetName();
445       } 
446       Printf("%s",str.Data());       
447     }//if(tid==-1)      
448   }//events loop
449   pAL->UnloadHeader();  pAL->UnloadKinematics();
450   return mtid;
451 }
452 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
453 Int_t AliRICHParam::StackCount(Int_t pid,Int_t evt)
454 {
455 // Counts total number of particles of given sort (including secondary) for a given event
456   AliRunLoader *pAL=AliRunLoader::Open(); 
457   pAL->GetEvent(evt);    
458   if(pAL->LoadHeader()) return 0;
459   if(pAL->LoadKinematics()) return 0;
460   AliStack *pStack=pAL->Stack();
461   
462   Int_t iCnt=0;
463   for(Int_t i=0;i<pStack->GetNtrack();i++) if(pStack->Particle(i)->GetPdgCode()==pid) iCnt++;
464   
465   pAL->UnloadHeader();  pAL->UnloadKinematics();
466   return iCnt;
467 }
468 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++