]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHParam.cxx
Check the Minuit status after each fit
[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 Bool_t         AliRICHParam::fgIsWireSag            =kTRUE;   //take ware sagita into account?
40 Bool_t         AliRICHParam::fgIsResolveClusters    =kTRUE;   //do cluster resolving?
41 Bool_t         AliRICHParam::fgIsFeedback           =kTRUE;   //generate feedback photons?
42 Bool_t         AliRICHParam::fgIsTestBeam           =kFALSE;  //special test beam configuration
43
44 Int_t    AliRICHParam::fgHV[kNsectors]        ={2050,2050,2050,2050,2050,2050};
45 Int_t    AliRICHParam::fgNsigmaTh             =4;
46 Float_t  AliRICHParam::fgSigmaThMean          =1.132; //QDC 
47 Float_t  AliRICHParam::fgSigmaThSpread        =0.035; //     
48 Double_t AliRICHParam::fgErrChrom[4][330];                       //
49 Double_t AliRICHParam::fgErrGeom[4][330];                        //
50 Double_t AliRICHParam::fgErrLoc[4][330];                         //Chromatic, Geometric and Localization array to parametrize SigmaCerenkov
51 Double_t AliRICHParam::fgMass[5]              ={0.00051,0.10566,0.13957,0.49360,0.93828};  
52
53
54 Double_t AliRICHParam::fEckovMin=5.5e-9; //GeV
55 Double_t AliRICHParam::fEckovMax=8.5e-9; //GeV
56 TF1 AliRICHParam::fgAbsC6F14("RabsC4F14","6512.39*(x<=7.75e-9)+(x>7.75e-9)*0.039/(-0.166+0.063e9*x-8.01e7*x^2+3.39e5*x^3)" ,fEckovMin,fEckovMax); 
57 TF1 AliRICHParam::fgAbsSiO2 ("RabsSiO2" ,"333"                                                                             ,fEckovMin,fEckovMax); 
58 TF1 AliRICHParam::fgAbsCH4  ("RabsCH4"  ,"6512.39*(x<=7.75e-9)+(x>7.75e-9)*0.039/(-0.166+0.063e9*x-8.01e7*x^2+3.39e5*x^3)" ,fEckovMin,fEckovMax);
59 TF1 AliRICHParam::fgAbsAir  ("RabsAir"  ,"500"                                                                             ,fEckovMin,fEckovMax);  //len ???
60
61 TF1 AliRICHParam::fgIdxAir  ("RidxAir"  ,"1+1e-8*(8342.13 + 2406030/(130-(1.23984e-9/x)^2)+15597/(38.9-(1.23984e-9/x)^2))" ,fEckovMin,fEckovMax);  //???
62 TF1 AliRICHParam::fgIdxSiO2 ("RidxSiO2" ,"sqrt(1+46.411/(10.666*10.666-x*x*1e18)+228.71/(18.125*18.125-x*x*1e18))"         ,fEckovMin,fEckovMax);  //TDR p.35
63 TF1 AliRICHParam::fgIdxCH4  ("RidxCH4"  ,"1+0.12489e-6/(2.62e-4 - (1239.84e-9/x)^-2)"                                      ,fEckovMin,fEckovMax);  //Olav preprint
64 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
65 AliRICHParam::AliRICHParam():TNamed("RichParam","default version") 
66 {
67 // Here all the intitializition is taken place when AliRICHParam::Instance() is invoked for the first time.
68 // In particulare, matrices to be used for LORS<->MARS trasnformations are initialized from TGeo structure.    
69 // Note that TGeoManager should be already initialized from geometry.root file  
70   for(Int_t iCh=0;iCh<kNchambers;iCh++) fMatrix[iCh]=(TGeoHMatrix*)gGeoManager->GetVolume("ALIC")->GetNode(Form("RICH_%i",iCh+1))->GetMatrix();
71   CdbRead(0,0);
72   fgInstance=this; 
73 }//ctor
74 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
75 Float_t AliRICHParam::AbsCH4(Float_t eV)
76 {
77 // Evaluate the absorbtion lenght of CH4 for a photon of energy eV in electron-volts
78   const Float_t kLoschmidt=2.686763e19;                                      // LOSCHMIDT NUMBER IN CM-3
79   const Float_t kPressure=750.0;          //mm of Hg
80   const Float_t kTemperature=283.0;       //K (10 grad C)                               
81   const Float_t kPn=kPressure/760.;
82   const Float_t kTn=kTemperature/273.16;
83   const Float_t kC0=-1.655279e-1;
84   const Float_t kC1= 6.307392e-2;
85   const Float_t kC2=-8.011441e-3;
86   const Float_t kC3= 3.392126e-4;
87                 
88   Float_t crossSection=0;                        
89   if (eV<7.75) 
90     crossSection=0.06e-22;
91   else                 //------ METHANE CROSS SECTION cm-2 ASTROPH. J. 214, L47 (1978)                                               
92     crossSection=(kC0+kC1*eV+kC2*eV*eV+kC3*eV*eV*eV)*1.e-18;
93     
94     Float_t density=kLoschmidt*kPn/kTn; //CH4 molecular concentration (cm^-3)
95     return 1.0/(density*crossSection);
96 }//AbsoCH4()
97 //__________________________________________________________________________________________________sss
98 void AliRICHParam::CdbRead(Int_t run,Int_t version)
99 {
100 // This methode read all the calibration information and initialise corresponding fields for requested run number
101 // Arguments: run - run number for which to retrieve calibration
102 //            version- version number   
103 //   Returns: none      
104
105   AliCDBEntry *pEntry=AliCDBManager::Instance()->Get("RICH/RICHConfig/RefIdxC6F14",run,0,version); //try to get from common local storage  
106   if(pEntry){
107     fIdxC6F14=(TF2*)pEntry->GetObject(); delete pEntry;
108   }else{
109     AliWarning("No valid calibarion, the hardcoded will be used!");
110     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
111     fIdxC6F14->SetUniqueID(20);//T=20 deg C
112   }
113 }//CdbRead()
114 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
115 void AliRICHParam::Print(Option_t* opt) const
116 {
117 // print some usefull (hopefully) info on some internal guts of RICH parametrisation 
118   Printf("Pads in chamber (%3i,%3i) in sector (%2i,%2i) pad size (%4.2f,%4.2f)",NpadsX(),NpadsY(),NpadsXsec(),NpadsYsec(),PadSizeX(),PadSizeY());
119   Printf("Resolve clusters %i sagita %i",IsResolveClusters(),IsWireSag()); 
120   
121   for(Int_t i=0;i<kNchambers;i++) fMatrix[i]->Print(opt);
122 }//Print()
123 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
124 void AliRICHParam::TestSeg()
125 {
126 // Provides a set of pictures to test segementation currently in use.    
127 // Arguments: none
128 //   Returns: none    
129   new TCanvas("pads","PC segmentation - pads display",700,600);
130   gPad->Range(-5,-5,PcSizeX()+5,PcSizeY()+15);
131   TVector p(2);   TVector2 c;    TVector2 b;   //current: pad, pad center, pad boundary
132 // list of corners:
133   Double_t x0=0,x1=SecSizeX(),x2=SecSizeX()+DeadZone()                                                      ,x3=PcSizeX();  
134   Double_t y0=0,y1=SecSizeY(),y2=SecSizeY()+DeadZone(),y3=2*SecSizeY()+DeadZone(),y4=PcSizeY()-SecSizeY(),y5=PcSizeY();  
135   DrawSectors();
136 //header 
137   TLatex t;
138   t.SetTextSize(0.02); t.SetTextColor(kBlack); t.SetTextAlign(11);
139   t.DrawLatex(0,PcSizeY()+10,Form("IP in front of this page. pad size %.2fx%.2fcm   dead zone %.2fcm",PadSizeX(),PadSizeY(),DeadZone()));
140   t.DrawLatex(0,PcSizeY()+ 5,Form("Pc  %.2fx%.2f cm %ix%i pads               Sec %.2fx%.2f cm %ix%i pads",
141                                          PcSizeX()     , PcSizeY()     , NpadsX()    , NpadsY()                                 ,
142                                          SecSizeX() , SecSizeY() , NpadsXsec() , NpadsYsec()                              ));
143 //sectors  
144   t.SetTextSize(0.015); t.SetTextColor(kRed); t.SetTextAlign(22);
145   c=Pad2Loc( 40, 24); t.DrawText(c.X(),c.Y(),Form("sec 1 (%.2f,%.2f)",c.X(),c.Y()  ));  
146   c=Pad2Loc( 40, 75); t.DrawText(c.X(),c.Y(),Form("sec 3 (%.2f,%.2f)",c.X(),c.Y()  ));  
147   c=Pad2Loc( 40,121); t.DrawText(c.X(),c.Y(),Form("sec 5 (%.2f,%.2f)",c.X(),c.Y()  ));  
148   c=Pad2Loc(120, 24); t.DrawText(c.X(),c.Y(),Form("sec 2 (%.2f,%.2f)",c.X(),c.Y()  ));  
149   c=Pad2Loc(120, 75); t.DrawText(c.X(),c.Y(),Form("sec 4 (%.2f,%.2f)",c.X(),c.Y()  ));  
150   c=Pad2Loc(120,121); t.DrawText(c.X(),c.Y(),Form("sec 6 (%.2f,%.2f)",c.X(),c.Y()  ));  
151 //coners  
152   t.SetTextSize(0.015); t.SetTextColor(kBlue);
153
154   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()));
155   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()));
156   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()));
157   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()));
158   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()));
159   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()));
160   
161   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()));
162   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()));
163   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()));
164   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()));
165   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()));
166   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()));
167   
168   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()));
169   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()));
170   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()));
171   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()));
172   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()));
173   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()));
174   
175   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()));
176   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()));
177   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()));
178   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()));
179   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()));
180   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()));
181 }//TestSeg()
182 //__________________________________________________________________________________________________
183 void AliRICHParam::TestResp()
184 {
185 // Provides a set of plot to check the response parametrisation currently in use.  
186 // Arguments: none
187 //   Returns: none    
188   TCanvas *pC=new TCanvas("c","Amplification test",900,800);
189   pC->Divide(1,2);
190   
191   
192   const Int_t kNpoints=8;
193   THStack *pStackPhot=new THStack("StackPhot","photons");
194   THStack *pStackMip =new THStack("StackMip","mips");
195   TLegend *pLeg=new TLegend(0.6,0.2,0.9,0.5,"legend");    
196   TH1F *apHphot[kNpoints];
197   TH1F *apHmip[kNpoints];
198   
199   Double_t starty=0;
200   Double_t deltay=AliRICHParam::SecSizeY()/kNpoints;
201   
202   for(int i=0;i<kNpoints;i++){
203     apHphot[i]=new TH1F(Form("hphot%i",i),"Qdc for Photon;QDC;Counts",500,0,500); apHphot[i]->SetLineColor(i);pStackPhot->Add(apHphot[i]);                 
204     apHmip[i] =new TH1F(Form("hmip%i",i),"Qdc for Mip;QDC;Counts",4000,0,4000);   apHmip[i]->SetLineColor(i);pStackMip->Add(apHmip[i]);                 
205     
206     pLeg->AddEntry(apHphot[i],Form("@(10,%5.2f->%5.2f)",starty+i*deltay,starty+i*deltay-SecSizeY()/2));
207   }
208         
209   
210   TVector2 x2(0,0);  
211   for(Int_t i=0;i<10000;i++){//events loop
212     for(int j=0;j<kNpoints;j++){
213       x2.Set(10,starty+j*deltay);
214       apHphot[j]->Fill(TotQdc(x2,0));
215       apHmip[j]->Fill(TotQdc(x2,gRandom->Landau(600,150)*1e-9));
216     }
217   }
218   
219   pC->cd(1);  pStackMip->Draw("nostack");
220   pC->cd(2);  pStackPhot->Draw("nostack"); pLeg->Draw();
221 }//TestResp()
222 //__________________________________________________________________________________________________
223 void AliRICHParam::TestTrans()
224 {
225 // Tests transformation methods
226 // Arguments: none
227 //   Returns: none    
228   
229   AliRICHParam *pParam=AliRICHParam::Instance();
230   Int_t iNpointsX=50,iNpointsY=50;  
231   new TCanvas("trasform","Test LORS-MARS transform"); TLatex t; t.SetTextSize(0.02);
232   
233   TView *pView=new TView(1);  pView->SetRange(-400,-400,-400,400,400,400);
234   DrawAxis();  
235   for(Int_t iCham=1;iCham<=7;iCham++){//chamber loop
236     AliRICHHelix helix(2.5,Norm(iCham).Theta()*TMath::RadToDeg(),Norm(iCham).Phi()*TMath::RadToDeg());
237     helix.RichIntersect(AliRICHParam::Instance());
238     TPolyMarker3D *pChamber=new TPolyMarker3D(iNpointsX*iNpointsY);
239     Int_t i=0;
240     for(Double_t x=0;x<PcSizeX();x+=PcSizeX()/iNpointsX)
241       for(Double_t y=0;y<PcSizeY();y+=PcSizeY()/iNpointsY){//step loop
242         TVector3 v3=pParam->Lors2Mars(iCham,x,y,kPc); TVector2 v2=pParam->Mars2Lors(iCham,v3,kPc);//LORS->MARS->LORS
243         Double_t dx=v2.X()-x , dy=v2.Y()-y; 
244         if(dx>0.000001 || dy>0.000001) Printf("Problem in MARS<->LORS transformations dx=%f dy=%f!!!",dx,dy);
245         pChamber->SetPoint(i++,v3.X(),v3.Y(),v3.Z());//Pc plane point in MARS
246       }//step loop
247     pChamber->SetMarkerSize(1);
248     pChamber->SetMarkerColor(iCham);
249     pChamber->Draw();  
250     helix.Draw();
251     t.SetNDC();t.SetTextColor(iCham); t.DrawText(0.1,iCham*0.1,Form("Chamber %i",iCham));    
252   }//chambers loop   
253 }//TestTrans()
254 //__________________________________________________________________________________________________
255 void AliRICHParam::DrawAxis()
256 {
257 // This utility methode draws axis on geometry scene  
258 // Arguments: none
259 //   Returns: none    
260   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};  
261   TPolyLine3D *pXaxis=new TPolyLine3D(2,x);pXaxis->SetLineColor(kRed);   pXaxis->Draw();
262   TPolyLine3D *pYaxis=new TPolyLine3D(2,y);pYaxis->SetLineColor(kGreen); pYaxis->Draw();
263   TPolyLine3D *pZaxis=new TPolyLine3D(2,z);pZaxis->SetLineColor(kBlue);  pZaxis->Draw();  
264 }
265 //__________________________________________________________________________________________________
266 void AliRICHParam::DrawSectors() 
267
268 // Utility methode draws RICH chamber sectors on event display.
269 // Arguments: none
270 //   Returns: none      
271   Double_t xLeft[5]  = {0,0,SecSizeX(),SecSizeX(),0};
272   Double_t xRight[5] = {SecSizeX()+DeadZone(),SecSizeX()+DeadZone(),PcSizeX(),PcSizeX(),SecSizeX()+DeadZone()};
273   
274   Double_t yDown[5]   = {0,SecSizeY(),SecSizeY(),0,0};
275   Double_t yCenter[5] = {  SecSizeY()+DeadZone(),2*SecSizeY()+DeadZone(),2*SecSizeY()+DeadZone(),
276                            SecSizeY()+DeadZone(),SecSizeY()+DeadZone()};  
277   Double_t yUp[5]     = {2*SecSizeY()+2*DeadZone(),PcSizeY(),PcSizeY(),2*SecSizeY()+2*DeadZone(),2*SecSizeY()+2*DeadZone()};
278   
279   TPolyLine *sec1 = new TPolyLine(5,xLeft ,yDown);    sec1->SetLineColor(21);  sec1->Draw();
280   TPolyLine *sec2 = new TPolyLine(5,xRight,yDown);    sec2->SetLineColor(21);  sec2->Draw();
281   TPolyLine *sec3 = new TPolyLine(5,xLeft ,yCenter);  sec3->SetLineColor(21);  sec3->Draw();
282   TPolyLine *sec4 = new TPolyLine(5,xRight,yCenter);  sec4->SetLineColor(21);  sec4->Draw();
283   TPolyLine *sec5 = new TPolyLine(5,xLeft, yUp);      sec5->SetLineColor(21);  sec5->Draw();
284   TPolyLine *sec6 = new TPolyLine(5,xRight,yUp);      sec6->SetLineColor(21);  sec6->Draw();
285 }//DrawSectors()
286 //__________________________________________________________________________________________________
287 void AliRICHParam::ReadErrFiles()
288 {
289 // Read the three files corresponding to Chrom,Geom and Loc They are parameters of a polynomial of 6th order...  ????????? go to CDB?
290 // Arguments: none
291 //   Returns: none    
292   
293   static Bool_t count = kFALSE;
294   
295   Float_t c0,c1,c2,c3,c;
296   Float_t g0,g1,g2,g3,g;
297   Float_t l0,l1,l2,l3,l;
298   
299   FILE *pChromErr, *pGeomErr, *pLocErr;  
300
301   if(!count) {
302      AliInfoGeneral("ReadErrFiles","reading RICH error parameters...");
303      pChromErr = fopen(Form("%s/RICH/RICHConfig/SigmaChromErr.txt",gSystem->Getenv("ALICE_ROOT")),"r");
304      pGeomErr  = fopen(Form("%s/RICH/RICHConfig/SigmaGeomErr.txt",gSystem->Getenv("ALICE_ROOT")),"r");
305      pLocErr   = fopen(Form("%s/RICH/RICHConfig/SigmaLocErr.txt",gSystem->Getenv("ALICE_ROOT")),"r");
306      if(!pChromErr||!pGeomErr||!pLocErr) {AliErrorGeneral("ReadErrFiles"," RICH ERROR READING Parameter FILES: can't open files!!! ");return;}
307      for(Int_t i=0;i<330;i++) {
308        fscanf(pChromErr,"%f%f%f%f%f\n",&c0,&c1,&c2,&c3,&c);
309        fscanf(pGeomErr,"%f%f%f%f%f\n",&g0,&g1,&g2,&g3,&g);
310        fscanf(pLocErr,"%f%f%f%f%f\n",&l0,&l1,&l2,&l3,&l);
311        fgErrChrom[0][i] = c0;
312        fgErrChrom[1][i] = c1;
313        fgErrChrom[2][i] = c2;
314        fgErrChrom[3][i] = c3;   
315        fgErrGeom[0][i] = g0;
316        fgErrGeom[1][i] = g1;
317        fgErrGeom[2][i] = g2;
318        fgErrGeom[3][i] = g3;    
319        fgErrLoc[0][i] = l0;
320        fgErrLoc[1][i] = l1;
321        fgErrLoc[2][i] = l2;
322        fgErrLoc[3][i] = l3;     
323      }
324      AliInfoGeneral("ReadErrFiles","DONE successfully!");
325      fclose(pChromErr);
326      fclose(pGeomErr);
327      fclose(pLocErr);
328   }
329   count = kTRUE;
330 }//ReadErrFiles()
331 //__________________________________________________________________________________________________
332 TVector3 AliRICHParam::SigmaSinglePhoton(Int_t partID, Double_t mom, Double_t theta, Double_t phi)
333
334 {
335 // Find sigma for single photon. It returns the thrree different errors. If you want
336 // to have the error---> TVector3.Mag()
337 // partID = 0,1,2,3,4 ---> e,mu,pi,k,p in agreement with AliPID
338   TVector3 v(-999,-999,-999);
339   Double_t pmom;
340
341   ReadErrFiles();
342   Double_t mass = fgMass[partID];
343   Double_t massRef = fgMass[4]; // all the files are calculated for protons...so mass ref is proton mass
344   pmom = mom*massRef/mass; // normalized momentum respect to proton...
345   if(pmom>PmodMax()) pmom = PmodMax();
346   Double_t oneOverRefIndex = 1/IdxC6F14(EckovMean());
347   Double_t pmin = mass*oneOverRefIndex/TMath::Sqrt(1-oneOverRefIndex*oneOverRefIndex);
348   if(pmom<pmin) return v;
349   Double_t Theta = theta*TMath::RadToDeg();
350   if(phi<0) phi+=TMath::TwoPi();
351   Double_t Phi = phi*TMath::RadToDeg();
352   v.SetX(Interpolate(fgErrChrom,pmom,Theta,Phi));
353   v.SetY(Interpolate(fgErrGeom,pmom,Theta,Phi));
354   v.SetZ(Interpolate(fgErrLoc,pmom,Theta,Phi));
355 //  v*=1.5; // take into account bigger errors due to multiplicity...to change in future
356
357   return v;
358 }//SigmaSinglePhoton
359 //__________________________________________________________________________________________________
360 TVector3 AliRICHParam::SigmaSinglePhoton(Double_t thetaCer, Double_t theta, Double_t phi)
361
362 {
363 // Find sigma for single photon. It returns the thrree different errors. If you want
364 // to have the error---> TVector3.Mag()
365 // partID = 0,1,2,3,4 ---> e,mu,pi,k,p in agreement with AliPID
366   TVector3 v(-999,-999,-999);
367   Double_t pmom;
368
369   ReadErrFiles();
370   Double_t massRef = fgMass[4]; // all the files are calculated for protons...so mass ref is proton mass
371   Double_t beta=1./(IdxC6F14(EckovMean())*TMath::Cos(thetaCer));
372   if(beta>=1) {
373     pmom=6.5; // above physical limi the error is calculated at the saturation...
374   } else {
375     Double_t gamma=1./TMath::Sqrt(1-beta*beta);
376     pmom = beta*gamma*massRef; // normalized momentum respect to proton...
377   }
378   if(pmom>PmodMax()) pmom = PmodMax();
379   Double_t oneOverRefIndex = 1/IdxC6F14(EckovMean());
380   Double_t pmin = massRef*oneOverRefIndex/TMath::Sqrt(1-oneOverRefIndex*oneOverRefIndex);
381   if(pmom<pmin) return v;
382   Double_t Theta = theta*TMath::RadToDeg();
383   if(phi<0) phi+=TMath::TwoPi();
384   Double_t Phi = phi*TMath::RadToDeg();
385   v.SetX(Interpolate(fgErrChrom,pmom,Theta,Phi));
386   v.SetY(Interpolate(fgErrGeom,pmom,Theta,Phi));
387   v.SetZ(Interpolate(fgErrLoc,pmom,Theta,Phi));
388 //  v*=1.5; // take into account bigger errors due to multiplicity...to change in future
389
390   return v;
391 }//SigmaSinglePhoton
392 //__________________________________________________________________________________________________
393 Double_t AliRICHParam::Interpolate(Double_t par[4][330], Double_t x, Double_t y, Double_t phi)
394   
395 {
396   static Double_t amin = 1.15; static Double_t astep  = 0.2;
397   static Double_t bmin = 0; static Double_t bstep  = 1;
398   
399   Double_t Phi = (phi - 180)/300.;
400   
401   Double_t Sigma[30][11];
402   
403   for(Int_t j=0;j<11;j++) { for(Int_t i=0;i<30;i++) {
404     Sigma[i][j] = par[0][j+11*i] + par[1][j+11*i]*Phi*Phi + par[2][j+11*i]*TMath::Power(Phi,4) + par[3][j+11*i]*TMath::Power(Phi,6);
405     }
406   }
407
408   Int_t i=0;Int_t j=0;
409   
410   i = (Int_t)((x-amin)/astep);
411   j = (Int_t)((y-bmin)/bstep);
412   Double_t ai = amin+i*astep;
413   Double_t ai1 = ai+astep;
414   Double_t bj = bmin+j*bstep;
415   Double_t bj1 = bj+bstep;
416   Double_t t = (x-ai)/(ai1-ai);
417   Double_t gj = (1-t)*Sigma[i][j]+t*Sigma[i+1][j];
418   Double_t gj1 = (1-t)*Sigma[i][j+1]+t*Sigma[i+1][j+1];
419   Double_t u = (y-bj)/(bj1-bj);
420   return (1-u)*gj+u*gj1;
421 }//Interpolate
422 //__________________________________________________________________________________________________
423 TVector3 AliRICHParam::ForwardTracing(TVector3 entranceTrackPoint, TVector3 vectorTrack, Double_t thetaC, Double_t phiC)
424 {
425 // Trace a single Ckov photon from a given emission point up to photocathode taking into account ref indexes of materials it travereses
426   TVector3 vBad(-999,-999,-999);
427   TVector3 nPlane(0,0,1);
428   Double_t planeZposition = 0.5*RadThick();
429   TVector3 planePoint(0,0,0.5*RadThick()); //this is plane parallel to window which contains emission point 
430   TVector3 emissionPoint = PlaneIntersect(vectorTrack,entranceTrackPoint,nPlane,planePoint);
431   Double_t thetaout,phiout;
432   AnglesInDRS(vectorTrack.Theta(),vectorTrack.Phi(),thetaC,phiC,thetaout,phiout);
433   TVector3 vectorPhotonInC6F14;  
434   vectorPhotonInC6F14.SetMagThetaPhi(1,thetaout,phiout);
435   planeZposition=RadThick();
436   planePoint.SetXYZ(0,0,planeZposition);
437   TVector3 entranceToSiO2Point = PlaneIntersect(vectorPhotonInC6F14,emissionPoint,nPlane,planePoint);
438
439   Double_t photonEn = EckovMean();
440   Double_t angleInSiO2 = SnellAngle(IdxC6F14(EckovMean()),IdxSiO2(EckovMean()),vectorPhotonInC6F14.Theta());if(angleInSiO2<0) return vBad;
441   TVector3 vectorPhotonInSiO2;
442   vectorPhotonInSiO2.SetMagThetaPhi(1,angleInSiO2,phiout);
443 //  planeZposition+=AliRICHParam::SiO2Thickness();
444   planeZposition+=WinThick();
445   planePoint.SetXYZ(0,0,planeZposition);
446   TVector3 entranceToCH4 = PlaneIntersect(vectorPhotonInSiO2,entranceToSiO2Point,nPlane,planePoint);
447 //  entranceToCH4.Dump();
448
449   //  Double_t angleInCH4 = SnellAngle(AliRICHParam::IndOfRefSiO2(6.755),AliRICHParam::IndOfRefCH4,angleInSiO2);
450   Double_t angleInCH4 = SnellAngle(IdxSiO2(photonEn),IdxCH4(photonEn),vectorPhotonInSiO2.Theta());if(angleInCH4<0) return vBad;
451   TVector3 vectorPhotonInCH4;
452   vectorPhotonInCH4.SetMagThetaPhi(1,angleInCH4,phiout);
453 //  planeZposition+=AliRICHParam::GapProx();
454   planeZposition+=Pc2Win();
455   planePoint.SetXYZ(0,0,planeZposition);
456   TVector3 impactToPC = PlaneIntersect(vectorPhotonInCH4,entranceToCH4,nPlane,planePoint);
457 //  impactToPC.Dump();
458   return impactToPC;
459 }//FowardTracing
460 //__________________________________________________________________________________________________
461 TVector3 AliRICHParam::PlaneIntersect(const TVector3 &lineDir,const TVector3 &linePoint,const TVector3 &planeNorm,const TVector3 &planePoint)
462 {
463 // Finds an intersection point between a line and plane.
464 // Arguments:  lineDir,linePoint    - vector along the line and any point of the line
465 //             planeNorm,planePoint - vector normal to the plane and any point of the plane
466 //   Returns:  point of intersection if any
467   if(planeNorm*lineDir==0) return   TVector3(-999,-999,-999);
468   TVector3 diff=planePoint-linePoint;
469   Double_t sint=(planeNorm*diff)/(planeNorm*lineDir);
470   return linePoint+sint*lineDir;
471 }//PlaneIntersect
472 //__________________________________________________________________________________________________ 
473 Double_t AliRICHParam::SnellAngle(Float_t n1, Float_t n2, Float_t theta1)
474 {
475 // Compute the angle of refraction out of Snell law
476 // Arguments: n1 - ref idx of first substance  
477 //            n2 - ref idx of second substance  
478 //            n1 - photon impact angle in the first substance i.e. angle between the photon direction and vector normal to the surface (radians)
479 //   Returns: photon refraction angle, i.e. angle in the second substance (radians)
480   Double_t sinref=(n1/n2)*TMath::Sin(theta1);
481   if(sinref>1.)    return -999;  
482   else             return TMath::ASin(sinref);
483 }//SnellAngle
484 //__________________________________________________________________________________________________
485 void AliRICHParam::AnglesInDRS(Double_t trackTheta,Double_t trackPhi,Double_t thetaCerenkov,Double_t phiCerenkov,Double_t &tout,Double_t &pout)
486 {
487 // Setup the rotation matrix of the track...
488
489   TRotation mtheta;
490   TRotation mphi;
491   TRotation minv;
492   TRotation mrot;
493   
494   mtheta.RotateY(trackTheta);
495   mphi.RotateZ(trackPhi);
496   
497   mrot = mphi * mtheta;
498     //  minv = mrot.Inverse();
499
500   TVector3 photonInRadiator(1,1,1);
501
502   photonInRadiator.SetTheta(thetaCerenkov);
503   photonInRadiator.SetPhi(phiCerenkov);
504   photonInRadiator = mrot * photonInRadiator;
505   tout=photonInRadiator.Theta();
506   pout=photonInRadiator.Phi();
507 }//AnglesInDRS
508 /*
509 void DrawRing()
510 {
511
512   //  Float_t xGraph[1000],yGraph[1000];
513
514   Float_t type;
515   Float_t MassOfParticle;
516   Float_t beta;
517   Float_t nfreon;
518
519   Float_t ThetaCerenkov;
520
521   Float_t Xtoentr = GetEntranceX();
522   Float_t Ytoentr = GetEntranceY();
523
524   Float_t pmod = GetTrackMomentum();
525   Float_t TrackTheta = GetTrackTheta();
526   Float_t TrackPhi = GetTrackPhi();
527
528   SetPhotonEnergy(AliRICHParam::MeanCkovEnergy());
529   SetFreonRefractiveIndex();
530
531   SetEmissionPoint(RadiatorWidth/2.);
532
533   ThetaCerenkov = GetThetaCerenkov();
534   FindBetaFromTheta(ThetaCerenkov);
535   nfreon = GetFreonRefractiveIndex();
536   
537   Int_t nPoints = 100;
538
539   Int_t nPointsToDraw = 0;
540   for(Int_t i=0;i<nPoints;i++)
541     {
542       Float_t phpad = 2*TMath::Pi()*i/nPoints;
543       SetThetaPhotonInTRS(thetacer);
544       SetPhiPhotonInTRS(phpad);
545       FindPhotonAnglesInDRS();
546       Float_t Radius = FromEmissionToCathode();
547       if (Radius == 999.) continue;
548       xGraph[nPointsToDraw] = GetXPointOnCathode() + GetShiftX();
549       yGraph[nPointsToDraw] = GetYPointOnCathode() + GetShiftY();
550       nPointsToDraw++;
551     }
552   gra = new TGraph(nPointsToDraw,xGraph,yGraph);
553   gra->Draw("AC"); 
554 }
555 //__________________________________________________________________________________________________
556 */
557 void AliRICHParam::TestHit2SDigs(Double_t x,Double_t y,Double_t e,Bool_t isNew)
558 {
559 //Test  hit->sdigits procedures
560 //Arguments: isNew - if true use new (abs pad) procedure else use old one (TVector)
561 //  Returns: none
562   TClonesArray *pSDigLst=new TClonesArray("AliRICHDigit");
563   Int_t iQtot=-1;
564   if(isNew){
565     iQtot=Hit2SDigs(10101,e,pSDigLst);        //new technique
566   }else{
567     iQtot=Hit2SDigs(TVector2(x,y),e,pSDigLst);//old technique
568   }
569   pSDigLst->Print();
570   Double_t dQsum=0;
571   for(Int_t i=0;i<pSDigLst->GetEntriesFast();i++)
572     dQsum+=((AliRICHDigit*)pSDigLst->At(i))->Qdc();
573   Printf("Qtot=%i Qsum=%.2f ",iQtot,dQsum);
574 }
575 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
576 Int_t AliRICHParam::Stack(Int_t evt,Int_t tid)
577 {
578 // Prints some usefull info from stack
579 // Arguments: evt - event number. if not -1 print info only for that event
580 //            tid - track id. if not -1 then print it and all it's mothers if any   
581 //   Returns: mother tid of the given tid if any
582   AliRunLoader *pAL=AliRunLoader::Open(); 
583   if(pAL->LoadHeader()) return -1;
584   if(pAL->LoadKinematics()) return -1;
585   
586   Int_t mtid=-1;
587   Int_t iNevt=pAL->GetNumberOfEvents();     Printf("This session contains %i event(s)",iNevt);
588   
589   for(Int_t iEvt=0;iEvt<iNevt;iEvt++){//events loop
590     if(evt!=-1 && evt!=iEvt) continue; //in case one needs to print the requested event, ignore all others
591     pAL->GetEvent(iEvt);    
592     AliStack *pStack=pAL->Stack();  
593     if(tid==-1){                        //print all tids for this event
594       for(Int_t i=0;i<pStack->GetNtrack();i++) pStack->Particle(i)->Print();
595       Printf("totally %i tracks including %i primaries for event %i out of %i event(s)",pStack->GetNtrack(),pStack->GetNprimary(),iEvt,iNevt);
596     }else{                              //print only this tid and it;s mothers
597       if(tid<0 || tid>pStack->GetNtrack()) {Printf("Wrong tid, valid tid range for event %i is 0-%i",iEvt,pStack->GetNtrack());break;}
598       TParticle *pTrack=pStack->Particle(tid); mtid=pTrack->GetFirstMother();
599       TString str=pTrack->GetName();
600       while((tid=pTrack->GetFirstMother()) >= 0){
601         pTrack=pStack->Particle(tid);
602         str+=" from ";str+=pTrack->GetName();
603       } 
604       Printf("%s",str.Data());       
605     }//if(tid==-1)      
606   }//events loop
607   pAL->UnloadHeader();  pAL->UnloadKinematics();
608   return mtid;
609 }
610 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
611 Int_t AliRICHParam::StackCount(Int_t pid,Int_t evt)
612 {
613 // Counts total number of particles of given sort (including secondary) for a given event
614   AliRunLoader *pAL=AliRunLoader::Open(); 
615   pAL->GetEvent(evt);    
616   if(pAL->LoadHeader()) return 0;
617   if(pAL->LoadKinematics()) return 0;
618   AliStack *pStack=pAL->Stack();
619   
620   Int_t iCnt=0;
621   for(Int_t i=0;i<pStack->GetNtrack();i++) if(pStack->Particle(i)->GetPdgCode()==pid) iCnt++;
622   
623   pAL->UnloadHeader();  pAL->UnloadKinematics();
624   return iCnt;
625 }
626 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++