// ************************************************************************** // * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * // * * // * Author: The ALICE Off-line Project. * // * Contributors are mentioned in the code where appropriate. * // * * // * Permission to use, copy, modify and distribute this software and its * // * documentation strictly for non-commercial purposes is hereby granted * // * without fee, provided that the above copyright notice appears in all * // * copies and that both the copyright notice and this permission notice * // * appear in the supporting documentation. The authors make no claims * // * about the suitability of this software for any purpose. It is * // * provided "as is" without express or implied warranty. * // ************************************************************************** #include "AliRICHParam.h" //class header #include "AliESD.h" #include //TestXXX() #include #include #include #include #include #include #include #include #include #include #include #include //CdbRead() #include //CdbRead() #include //CdbRead() #include //Stack() #include //Stack() #include //Stack() #include "AliRICHHelix.h" //TestTrans() #include ClassImp(AliRICHParam) AliRICHParam * AliRICHParam::fgInstance =0x0; //singleton pointer //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ AliRICHParam::AliRICHParam():TNamed("RichParam","default version") { // Here all the intitializition is taken place when AliRICHParam::Instance() is invoked for the first time. // In particulare, matrices to be used for LORS<->MARS trasnformations are initialized from TGeo structure. // Note that TGeoManager should be already initialized from geometry.root file for(Int_t iCh=0;iChGetVolume("ALIC")->GetNode(Form("RICH_%i",iCh+1))->GetMatrix(); CdbRead(0,0); fgInstance=this; }//ctor //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ Float_t AliRICHParam::AbsCH4(Float_t eV) { // Evaluate the absorbtion lenght of CH4 for a photon of energy eV in electron-volts const Float_t kLoschmidt=2.686763e19; // LOSCHMIDT NUMBER IN CM-3 const Float_t kPressure=750.0; //mm of Hg const Float_t kTemperature=283.0; //K (10 grad C) const Float_t kPn=kPressure/760.; const Float_t kTn=kTemperature/273.16; const Float_t kC0=-1.655279e-1; const Float_t kC1= 6.307392e-2; const Float_t kC2=-8.011441e-3; const Float_t kC3= 3.392126e-4; Float_t crossSection=0; if (eV<7.75) crossSection=0.06e-22; else //------ METHANE CROSS SECTION cm-2 ASTROPH. J. 214, L47 (1978) crossSection=(kC0+kC1*eV+kC2*eV*eV+kC3*eV*eV*eV)*1.e-18; Float_t density=kLoschmidt*kPn/kTn; //CH4 molecular concentration (cm^-3) return 1.0/(density*crossSection); }//AbsoCH4() //__________________________________________________________________________________________________sss void AliRICHParam::CdbRead(Int_t run,Int_t version) { // This methode read all the calibration information and initialise corresponding fields for requested run number // Arguments: run - run number for which to retrieve calibration // version- version number // Returns: none AliCDBEntry *pEntry=AliCDBManager::Instance()->Get("RICH/RICHConfig/RefIdxC6F14",run,0,version); //try to get from common local storage if(pEntry){ fIdxC6F14=(TF2*)pEntry->GetObject(); if(!(AliCDBManager::Instance()->GetCacheFlag())) delete pEntry; }else{ AliWarning("No valid calibarion, the hardcoded will be used!"); 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 fIdxC6F14->SetUniqueID(20);//T=20 deg C } }//CdbRead() //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ void AliRICHParam::Print(Option_t* opt) const { // print some usefull (hopefully) info on some internal guts of RICH parametrisation Printf("Pads in chamber (%3i,%3i) in sector (%2i,%2i) pad size (%4.2f,%4.2f)",NpadsX(),NpadsY(),NpadsXsec(),NpadsYsec(),PadSizeX(),PadSizeY()); for(Int_t i=0;iPrint(opt); }//Print() //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ void AliRICHParam::TestSeg() { // Provides a set of pictures to test segementation currently in use. // Arguments: none // Returns: none new TCanvas("pads","View from electronics side, IP is behind the picture."); gPad->Range(-5,-5,PcSizeX()+5,PcSizeY()+5); // list of corners: DrawSectors(kTRUE); TLatex t; t.SetTextSize(0.02); t.SetTextColor(kRed) ;t.SetTextAlign(11); t.DrawLatex(0 ,PcSizeY(),Form("Pad size %.2fx%.2fcm dead zone %.2fcm",PadSizeX(),PadSizeY(),DeadZone())); t.SetTextColor(kBlue) ;t.SetTextAlign(21); t.DrawLatex(SecSizeX(),PcSizeY(),Form("Pc size %.2fx%.2fcm %ix%i pads" ,PcSizeX() ,PcSizeY() ,NpadsX() ,NpadsY())); t.SetTextColor(kMagenta);t.SetTextAlign(31); t.DrawLatex(PcSizeX() ,PcSizeY(),Form("Sec size %.2fx%.2fcm %ix%i pads" ,SecSizeX(),SecSizeY(),NpadsXsec() ,NpadsYsec())); // TVector p(2); TVector2 c; TVector2 b; //current: pad, pad center, pad boundary // Double_t x0=0,x1=SecSizeX(),x2=SecSizeX()+DeadZone() ,x3=PcSizeX(); // Double_t y0=0,y1=SecSizeY(),y2=SecSizeY()+DeadZone(),y3=2*SecSizeY()+DeadZone(),y4=PcSizeY()-SecSizeY(),y5=PcSizeY(); t.SetTextSize(0.015); t.SetTextColor(kBlue); // 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())); // 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())); // 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())); // 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())); // 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())); // 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())); // // 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())); // 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())); // 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())); // 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())); // 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())); // 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())); // // 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())); // 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())); // 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())); // 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())); // 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())); // 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())); // // 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())); // 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())); // 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())); // 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())); // 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())); // 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())); //Now all chambers view TCanvas *pC=new TCanvas("cc","Chambers view from electronics side - IP is behind the picture"); pC->Divide(3,3); t.SetTextSize(0.05); for(Int_t i=1;i<=7;i++){ if(i==1) pC->cd(9); if(i==2) pC->cd(8); if(i==3) pC->cd(6); if(i==4) pC->cd(5); if(i==5) pC->cd(4); if(i==6) pC->cd(2); if(i==7) pC->cd(1); gPad->Range(-5,-5,PcSizeX()+5,PcSizeY()+10); DrawSectors(kTRUE); t.SetTextColor(kBlack); t.SetTextAlign(21); t.DrawText(SecSizeX(),PcSizeY(),Form("Cham %i",i)); t.SetTextColor(41); t.SetTextAlign(11); t.DrawText(0 ,PcSizeY(),Form("DDL ID %i",1536+2*i-2)); //left half of chamber t.SetTextColor(29); t.SetTextAlign(31); t.DrawText(PcSizeX() ,PcSizeY(),Form("DDL ID %i",1536+2*i-1)); //right half of chamber } }//TestSeg() //__________________________________________________________________________________________________ void AliRICHParam::TestResp() { // Provides a set of plot to check the response parametrisation currently in use. // Arguments: none // Returns: none TCanvas *pC=new TCanvas("c","Amplification test",900,800); pC->Divide(1,2); const Int_t kNpoints=8; THStack *pStackPhot=new THStack("StackPhot","photons"); THStack *pStackMip =new THStack("StackMip","mips"); TLegend *pLeg=new TLegend(0.6,0.2,0.9,0.5,"legend"); TH1F *apHphot[kNpoints]; TH1F *apHmip[kNpoints]; Double_t starty=0; Double_t deltay=AliRICHParam::SecSizeY()/kNpoints; for(int i=0;iSetLineColor(i);pStackPhot->Add(apHphot[i]); apHmip[i] =new TH1F(Form("hmip%i",i),"Qdc for Mip;QDC;Counts",4000,0,4000); apHmip[i]->SetLineColor(i);pStackMip->Add(apHmip[i]); pLeg->AddEntry(apHphot[i],Form("@(10,%5.2f->%5.2f)",starty+i*deltay,starty+i*deltay-SecSizeY()/2)); } TVector2 x2(0,0); for(Int_t i=0;i<10000;i++){//events loop for(int j=0;jFill(TotQdc(x2,0)); apHmip[j]->Fill(TotQdc(x2,gRandom->Landau(600,150)*1e-9)); } } pC->cd(1); pStackMip->Draw("nostack"); pC->cd(2); pStackPhot->Draw("nostack"); pLeg->Draw(); }//TestResp() //__________________________________________________________________________________________________ void AliRICHParam::TestTrans() { // Tests transformation methods // Arguments: none // Returns: none AliRICHParam *pParam=AliRICHParam::Instance(); Int_t iNpointsX=50,iNpointsY=50; new TCanvas("trasform","Test LORS-MARS transform"); TLatex t; t.SetTextSize(0.02); TView *pView=new TView(1); pView->SetRange(-400,-400,-400,400,400,400); DrawAxis(); for(Int_t iCham=1;iCham<=7;iCham++){//chamber loop AliRICHHelix helix(2.5,Norm(iCham).Theta()*TMath::RadToDeg(),Norm(iCham).Phi()*TMath::RadToDeg()); helix.RichIntersect(AliRICHParam::Instance()); TPolyMarker3D *pChamber=new TPolyMarker3D(iNpointsX*iNpointsY); Int_t i=0; for(Double_t x=0;xLors2Mars(iCham,x,y,kPc); TVector2 v2=pParam->Mars2Lors(iCham,v3,kPc);//LORS->MARS->LORS Double_t dx=v2.X()-x , dy=v2.Y()-y; if(dx>0.000001 || dy>0.000001) Printf("Problem in MARS<->LORS transformations dx=%f dy=%f!!!",dx,dy); pChamber->SetPoint(i++,v3.X(),v3.Y(),v3.Z());//Pc plane point in MARS }//step loop pChamber->SetMarkerSize(1); pChamber->SetMarkerColor(iCham); pChamber->Draw(); helix.Draw(); t.SetNDC();t.SetTextColor(iCham); t.DrawText(0.1,iCham*0.1,Form("Chamber %i",iCham)); }//chambers loop }//TestTrans() //__________________________________________________________________________________________________ void AliRICHParam::DrawAxis() { // This utility methode draws axis on geometry scene // Arguments: none // Returns: none 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}; TPolyLine3D *pXaxis=new TPolyLine3D(2,x);pXaxis->SetLineColor(kRed); pXaxis->Draw(); TPolyLine3D *pYaxis=new TPolyLine3D(2,y);pYaxis->SetLineColor(kGreen); pYaxis->Draw(); TPolyLine3D *pZaxis=new TPolyLine3D(2,z);pZaxis->SetLineColor(kBlue); pZaxis->Draw(); } //__________________________________________________________________________________________________ void AliRICHParam::DrawSectors(Bool_t isInfo) { // Utility methode draws RICH chamber sectors on event display. // Arguments: none // Returns: none Double_t xLeft[5] = {0,0,SecSizeX(),SecSizeX(),0}; Double_t xRight[5] = {SecSizeX()+DeadZone(),SecSizeX()+DeadZone(),PcSizeX(),PcSizeX(),SecSizeX()+DeadZone()}; Double_t yDown[5] = {0,SecSizeY(),SecSizeY(),0,0}; Double_t yCenter[5] = { SecSizeY()+DeadZone(),2*SecSizeY()+DeadZone(),2*SecSizeY()+DeadZone(), SecSizeY()+DeadZone(),SecSizeY()+DeadZone()}; Double_t yUp[5] = {2*SecSizeY()+2*DeadZone(),PcSizeY(),PcSizeY(),2*SecSizeY()+2*DeadZone(),2*SecSizeY()+2*DeadZone()}; Int_t iColorLeft=29,iColorRight=41; TPolyLine *sec[6]; sec[0]= new TPolyLine(5,xLeft ,yDown); sec[1]= new TPolyLine(5,xRight,yDown); sec[2]= new TPolyLine(5,xLeft ,yCenter); sec[3]= new TPolyLine(5,xRight,yCenter); sec[4]= new TPolyLine(5,xLeft, yUp); sec[5]= new TPolyLine(5,xRight,yUp); for(Int_t iSec=0;iSec<6;iSec++){ if(isInfo){ (iSec%2)? sec[iSec]->SetFillColor(iColorLeft): sec[iSec]->SetFillColor(iColorRight); sec[iSec]->Draw("f"); }else sec[iSec]->Draw(); } if(!isInfo) return; TText t; t.SetTextSize(0.02); t.DrawText(32,20,"Sec 1"); t.DrawText(98,20,"Sec 2"); t.DrawText(32,65,"Sec 3"); t.DrawText(98,65,"Sec 4"); t.DrawText(32,106,"Sec 5"); t.DrawText(98,106,"Sec 6"); Double_t x246=SecSizeX()+DeadZone(); Double_t y34=SecSizeY()+DeadZone(); Double_t y56=2*y34; for(Int_t iRow=1;iRow<=8;iRow++){//dilogic chip serves 8x6 pads Double_t y=(iRow-1)*6*PadSizeY(); TLine *pL1=new TLine(0 ,y ,SecSizeX(),y ); pL1->Draw(); t.SetTextAlign(31); t.DrawText(0 ,y ,Form("r%i",iRow) );//sec1 TLine *pL2=new TLine(x246,y ,PcSizeX() ,y ); pL2->Draw(); t.SetTextAlign(11); t.DrawText(PcSizeX(),y ,Form("r%i",25-iRow));//sec2 TLine *pL3=new TLine(0 ,y+y34,SecSizeX(),y+y34); pL3->Draw(); t.SetTextAlign(31); t.DrawText(0 ,y+y34,Form("r%i",iRow+8) );//sec3 TLine *pL4=new TLine(x246,y+y34,PcSizeX() ,y+y34); pL4->Draw(); t.SetTextAlign(11); t.DrawText(PcSizeX(),y+y34,Form("r%i",17-iRow));//sec4 TLine *pL5=new TLine(0 ,y+y56,SecSizeX(),y+y56); pL5->Draw(); t.SetTextAlign(31); t.DrawText(0 ,y+y56,Form("r%i",iRow+16));//sec5 TLine *pL6=new TLine(x246,y+y56,PcSizeX() ,y+y56); pL6->Draw(); t.SetTextAlign(11); t.DrawText(PcSizeX(),y+y56,Form("r%i",9-iRow) );//sec4 } for(Int_t iDilogic=1;iDilogic<=10;iDilogic++){ Double_t x=(iDilogic-1)*8*PadSizeX(); TLine *pL1=new TLine(x ,0,x ,SecSizeY()); pL1->Draw(); t.DrawText(x ,0,Form("d%i",11-iDilogic)); TLine *pL2=new TLine(x+x246,0,x+x246,SecSizeY()); pL2->Draw(); t.DrawText(x+x246,0,Form("d%i",11-iDilogic)); TLine *pL3=new TLine(x ,y34,x ,y34+SecSizeY()); pL3->Draw(); t.DrawText(x ,y34,Form("d%i",11-iDilogic)); TLine *pL4=new TLine(x+x246,y34,x+x246,y34+SecSizeY()); pL4->Draw(); t.DrawText(x+x246,y34,Form("d%i",11-iDilogic)); TLine *pL5=new TLine(x ,y56,x ,y56+SecSizeY()); pL5->Draw(); t.DrawText(x ,y56,Form("d%i",11-iDilogic)); TLine *pL6=new TLine(x+x246,y56,x+x246,y56+SecSizeY()); pL6->Draw(); t.DrawText(x+x246,y56,Form("d%i",11-iDilogic)); } t.SetTextAlign(13); t.DrawText(0 ,0,"pad1"); t.DrawText(x246 ,0,"pad81"); t.SetTextAlign(33); t.DrawText(SecSizeX(),0,"pad80"); t.DrawText(PcSizeX(),0,"pad160"); }//DrawSectors() //__________________________________________________________________________________________________ TVector3 AliRICHParam::ForwardTracing(TVector3 entranceTrackPoint, TVector3 vectorTrack, Double_t thetaC, Double_t phiC) { // Trace a single Ckov photon from a given emission point up to photocathode taking into account ref indexes of materials it travereses TVector3 vBad(-999,-999,-999); TVector3 nPlane(0,0,1); Double_t planeZposition = 0.5*RadThick(); TVector3 planePoint(0,0,0.5*RadThick()); //this is plane parallel to window which contains emission point TVector3 emissionPoint = PlaneIntersect(vectorTrack,entranceTrackPoint,nPlane,planePoint); Double_t thetaout,phiout; AnglesInDRS(vectorTrack.Theta(),vectorTrack.Phi(),thetaC,phiC,thetaout,phiout); TVector3 vectorPhotonInC6F14; vectorPhotonInC6F14.SetMagThetaPhi(1,thetaout,phiout); planeZposition=RadThick(); planePoint.SetXYZ(0,0,planeZposition); TVector3 entranceToSiO2Point = PlaneIntersect(vectorPhotonInC6F14,emissionPoint,nPlane,planePoint); Double_t photonEn = EckovMean(); Double_t angleInSiO2 = SnellAngle(IdxC6F14(EckovMean()),IdxSiO2(EckovMean()),vectorPhotonInC6F14.Theta());if(angleInSiO2<0) return vBad; TVector3 vectorPhotonInSiO2; vectorPhotonInSiO2.SetMagThetaPhi(1,angleInSiO2,phiout); // planeZposition+=AliRICHParam::SiO2Thickness(); planeZposition+=WinThick(); planePoint.SetXYZ(0,0,planeZposition); TVector3 entranceToCH4 = PlaneIntersect(vectorPhotonInSiO2,entranceToSiO2Point,nPlane,planePoint); // entranceToCH4.Dump(); // Double_t angleInCH4 = SnellAngle(AliRICHParam::IndOfRefSiO2(6.755),AliRICHParam::IndOfRefCH4,angleInSiO2); Double_t angleInCH4 = SnellAngle(IdxSiO2(photonEn),IdxCH4(photonEn),vectorPhotonInSiO2.Theta());if(angleInCH4<0) return vBad; TVector3 vectorPhotonInCH4; vectorPhotonInCH4.SetMagThetaPhi(1,angleInCH4,phiout); // planeZposition+=AliRICHParam::GapProx(); planeZposition+=Pc2Win(); planePoint.SetXYZ(0,0,planeZposition); TVector3 impactToPC = PlaneIntersect(vectorPhotonInCH4,entranceToCH4,nPlane,planePoint); // impactToPC.Dump(); return impactToPC; }//FowardTracing //__________________________________________________________________________________________________ TVector3 AliRICHParam::PlaneIntersect(const TVector3 &lineDir,const TVector3 &linePoint,const TVector3 &planeNorm,const TVector3 &planePoint) { // Finds an intersection point between a line and plane. // Arguments: lineDir,linePoint - vector along the line and any point of the line // planeNorm,planePoint - vector normal to the plane and any point of the plane // Returns: point of intersection if any if(planeNorm*lineDir==0) return TVector3(-999,-999,-999); TVector3 diff=planePoint-linePoint; Double_t sint=(planeNorm*diff)/(planeNorm*lineDir); return linePoint+sint*lineDir; }//PlaneIntersect //__________________________________________________________________________________________________ Double_t AliRICHParam::SnellAngle(Float_t n1, Float_t n2, Float_t theta1) { // Compute the angle of refraction out of Snell law // Arguments: n1 - ref idx of first substance // n2 - ref idx of second substance // n1 - photon impact angle in the first substance i.e. angle between the photon direction and vector normal to the surface (radians) // Returns: photon refraction angle, i.e. angle in the second substance (radians) Double_t sinref=(n1/n2)*TMath::Sin(theta1); if(sinref>1.) return -999; else return TMath::ASin(sinref); }//SnellAngle //__________________________________________________________________________________________________ void AliRICHParam::AnglesInDRS(Double_t trackTheta,Double_t trackPhi,Double_t thetaCerenkov,Double_t phiCerenkov,Double_t &tout,Double_t &pout) { // Setup the rotation matrix of the track... TRotation mtheta; TRotation mphi; TRotation minv; TRotation mrot; mtheta.RotateY(trackTheta); mphi.RotateZ(trackPhi); mrot = mphi * mtheta; // minv = mrot.Inverse(); TVector3 photonInRadiator(1,1,1); photonInRadiator.SetTheta(thetaCerenkov); photonInRadiator.SetPhi(phiCerenkov); photonInRadiator = mrot * photonInRadiator; tout=photonInRadiator.Theta(); pout=photonInRadiator.Phi(); }//AnglesInDRS /* void DrawRing() { // Float_t xGraph[1000],yGraph[1000]; Float_t type; Float_t MassOfParticle; Float_t beta; Float_t nfreon; Float_t ThetaCerenkov; Float_t Xtoentr = GetEntranceX(); Float_t Ytoentr = GetEntranceY(); Float_t pmod = GetTrackMomentum(); Float_t TrackTheta = GetTrackTheta(); Float_t TrackPhi = GetTrackPhi(); SetPhotonEnergy(AliRICHParam::MeanCkovEnergy()); SetFreonRefractiveIndex(); SetEmissionPoint(RadiatorWidth/2.); ThetaCerenkov = GetThetaCerenkov(); FindBetaFromTheta(ThetaCerenkov); nfreon = GetFreonRefractiveIndex(); Int_t nPoints = 100; Int_t nPointsToDraw = 0; for(Int_t i=0;iDraw("AC"); } //__________________________________________________________________________________________________ */ void AliRICHParam::TestHit2SDigs(Double_t x,Double_t y,Double_t e,Bool_t isNew) { //Test hit->sdigits procedures //Arguments: isNew - if true use new (abs pad) procedure else use old one (TVector) // Returns: none TClonesArray *pSDigLst=new TClonesArray("AliRICHDigit"); Int_t iQtot=-1; if(isNew){ iQtot=Hit2SDigs(10101,e,pSDigLst); //new technique }else{ iQtot=Hit2SDigs(TVector2(x,y),e,pSDigLst);//old technique } pSDigLst->Print(); Double_t dQsum=0; for(Int_t i=0;iGetEntriesFast();i++) dQsum+=((AliRICHDigit*)pSDigLst->At(i))->Qdc(); Printf("Qtot=%i Qsum=%.2f ",iQtot,dQsum); } //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ Int_t AliRICHParam::Stack(Int_t evt,Int_t tid) { // Prints some usefull info from stack // Arguments: evt - event number. if not -1 print info only for that event // tid - track id. if not -1 then print it and all it's mothers if any // Returns: mother tid of the given tid if any AliRunLoader *pAL=AliRunLoader::Open(); if(pAL->LoadHeader()) return -1; if(pAL->LoadKinematics()) return -1; Int_t mtid=-1; Int_t iNevt=pAL->GetNumberOfEvents(); Printf("This session contains %i event(s)",iNevt); for(Int_t iEvt=0;iEvtGetEvent(iEvt); AliStack *pStack=pAL->Stack(); if(tid==-1){ //print all tids for this event for(Int_t i=0;iGetNtrack();i++) pStack->Particle(i)->Print(); Printf("totally %i tracks including %i primaries for event %i out of %i event(s)",pStack->GetNtrack(),pStack->GetNprimary(),iEvt,iNevt); }else{ //print only this tid and it;s mothers if(tid<0 || tid>pStack->GetNtrack()) {Printf("Wrong tid, valid tid range for event %i is 0-%i",iEvt,pStack->GetNtrack());break;} TParticle *pTrack=pStack->Particle(tid); mtid=pTrack->GetFirstMother(); TString str=pTrack->GetName(); while((tid=pTrack->GetFirstMother()) >= 0){ pTrack=pStack->Particle(tid); str+=" from ";str+=pTrack->GetName(); } Printf("%s",str.Data()); }//if(tid==-1) }//events loop pAL->UnloadHeader(); pAL->UnloadKinematics(); return mtid; } //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ Int_t AliRICHParam::StackCount(Int_t pid,Int_t evt) { // Counts total number of particles of given sort (including secondary) for a given event AliRunLoader *pAL=AliRunLoader::Open(); pAL->GetEvent(evt); if(pAL->LoadHeader()) return 0; if(pAL->LoadKinematics()) return 0; AliStack *pStack=pAL->Stack(); Int_t iCnt=0; for(Int_t i=0;iGetNtrack();i++) if(pStack->Particle(i)->GetPdgCode()==pid) iCnt++; pAL->UnloadHeader(); pAL->UnloadKinematics(); return iCnt; } //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++