]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/Hdisp.C
unzoom set in display
[u/mrichter/AliRoot.git] / HMPID / Hdisp.C
1 TCanvas *pAll=0;
2 AliRunLoader *gAL=0; AliLoader *gHL=0; AliESD *gEsd=0; TTree *gEsdTr=0; AliHMPID *gH=0;
3 Int_t gEvt=0; Int_t gMaxEvt=0;
4 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
5 void Hdisp()
6 {//display events from files if any in current directory or simulated events
7   pAll=new TCanvas("all","",1000,900); pAll->Divide(3,3,0,0);pAll->ToggleEditor();
8   pAll->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",0,"","DoZoom(Int_t,Int_t,Int_t,TObject*)");
9   
10   if(gSystem->IsFileInIncludePath("galice.root")){// tries to open session
11     if(gAlice) delete gAlice;                                               //in case we execute this in aliroot delete default AliRun object 
12     gAL=AliRunLoader::Open();                                                                    //try to open galice.root from current dir 
13     gAL->LoadgAlice();                                                                           //take new AliRun object from galice.root   
14     gHL=gAL->GetDetectorLoader("HMPID");  gH=(AliHMPID*)gAL->GetAliRun()->GetDetector("HMPID");  //get HMPID object from galice.root
15     gMaxEvt=gAL->GetNumberOfEvents()-1;
16     gHL->LoadHits(); gHL->LoadDigits(); gHL->LoadRecPoints();
17
18     TFile *pEsdFl=TFile::Open("AliESDs.root"); gEsdTr=(TTree*) pEsdFl->Get("esdTree"); gEsdTr->SetBranchAddress("ESD", &gEsd);
19     pAll->cd(7); TButton *pBtn=new TButton("Next","ReadEvt()",0,0,0.2,0.1);   pBtn->Draw();
20     ReadEvt();
21   }else{
22     pAll->cd(7); TButton *pBtn=new TButton("Next","SimEvt()",0,0,0.2,0.1);   pBtn->Draw(); 
23     SimEvt();
24   }      
25 }
26 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 void ReadEvt()
28 {// Read curent event and display it assumes that session is alredy opened
29   TClonesArray hits("AliHMPIDHit");
30   if(gEvt>gMaxEvt) gEvt=0; if(gEvt<0) gEvt=gMaxEvt;
31   
32   gEsdTr->GetEntry(gEvt); gAL->GetEvent(gEvt); 
33   ReadHits(&hits); gHL->TreeD()->GetEntry(0); gHL->TreeR()->GetEntry(0);
34   
35   pAll->cd(3);  gPad->Clear(); TLatex txt;txt.DrawLatex(0.2,0.2,Form("Event %i (total %i)",gEvt,gMaxEvt));
36   DrawEvt(&hits,gH->DigLst(),gH->CluLst(),gEsd);
37   gEvt++;
38 }
39 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
40 void SimEvt()
41 {
42   TClonesArray hits("AliHMPIDHit"); 
43   TClonesArray sdig("AliHMPIDDigit"); 
44   TObjArray    digs(7); for(Int_t i=0;i<7;i++) digs.AddAt(new TClonesArray("AliHMPIDDigit"),i);
45   TObjArray    clus(7); for(Int_t i=0;i<7;i++) clus.AddAt(new TClonesArray("AliHMPIDCluster"),i);
46   AliESD esd;
47     
48   SimEsd(&esd);
49   SimHits(&esd,&hits);
50              AliHMPIDv1::Hit2Sdi(&hits,&sdig);                               
51       AliHMPIDDigitizer::Sdi2Dig(&sdig,&digs);     
52   AliHMPIDReconstructor::Dig2Clu(&digs,&clus);
53         AliHMPIDTracker::Recon(&esd,&clus);
54   
55   pAll->cd(3);  gPad->Clear(); TLatex txt;txt.DrawLatex(0.2,0.2,Form("Simulated event %i",gEvt));
56   DrawEvt(&hits,&digs,&clus,&esd);  
57   gEvt++;
58 }//SimEvt()
59 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
60 void SimEsd(AliESD *pEsd)
61 {
62   TParticle part; TLorentzVector mom;
63   for(Int_t iTrk=0;iTrk<25;iTrk++){//stack loop
64     part.SetPdgCode(kProton);
65     part.SetProductionVertex(0,0,0,0);  
66     Double_t eta= -0.4+gRandom->Rndm()*0.8; //rapidity is random [-0.4,+0.4]
67     Double_t phi= gRandom->Rndm()*1.4;      //phi is random      [ 0  , 80 ] degrees    
68     mom.SetPtEtaPhiM(2,eta,phi,part.GetMass());
69     part.SetMomentum(mom);
70     AliESDtrack trk(&part);
71     pEsd->AddTrack(&trk);
72   }//stack loop  
73 }//EsdFromStack()
74 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
75 void SimHits(AliESD *pEsd, TClonesArray *pHits)
76 {
77   AliHMPIDRecon rec;
78   const Int_t kCerenkov=50000050,kFeedback=50000051;
79   Int_t hc=0; TVector2 pos;
80   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop
81     AliESDtrack *pTrk=pEsd->GetTrack(iTrk);
82     Float_t xRa,yRa;
83     Int_t ch=AliHMPIDTracker::IntTrkCha(pTrk,xRa,yRa);
84     if(ch<0) continue; //this track does not hit HMPID
85     Float_t ckov=0.63;
86
87     Float_t th,ph,xPc,yPc,; pTrk->GetHMPIDtrk(xPc,yPc,th,ph); rec.SetTrack(xRa,yRa,th,ph); 
88     
89     if(!AliHMPIDDigit::IsInDead(xPc,yPc)) new((*pHits)[hc++]) AliHMPIDHit(ch,200e-9,kProton  ,iTrk,xPc,yPc);                 //mip hit
90     for(int i=0;i<4;i++)  new((*pHits)[hc++]) AliHMPIDHit(ch,7.5e-9,kFeedback,iTrk,gRandom->Rndm()*130,gRandom->Rndm()*126); //bkg hits 4 per track
91     for(int i=0;i<16;i++){
92       rec.TracePhot(ckov,gRandom->Rndm()*TMath::TwoPi(),pos);
93       new((*pHits)[hc++]) AliHMPIDHit(ch,7.5e-9,kCerenkov,iTrk,pos.X(),pos.Y());
94     }                      //photon hits  
95   }//tracks loop    
96 }//SimHits()
97 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
98 void ReadHits(TClonesArray *pHitLst)
99 {
100   pHitLst->Delete();  Int_t cnt=0;
101   for(Int_t iEnt=0;iEnt<gHL->TreeH()->GetEntries();iEnt++){       //TreeH loop
102     gHL->TreeH()->GetEntry(iEnt);                                 //get current entry (prim)                
103     for(Int_t iHit=0;iHit<gH->Hits()->GetEntries();iHit++){       //hits loop
104       AliHMPIDHit *pHit = (AliHMPIDHit*)gH->Hits()->At(iHit);     //get current hit        
105       new((*pHitLst)[cnt++]) AliHMPIDHit(*pHit);
106     }//hits loop for this entry
107   }//tree entries loop
108   
109 }//ReadHits()
110 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111 void DrawCh(Int_t iCh) 
112
113   gPad->Range(-10,-10,AliHMPIDDigit::SizeAllX()+5,AliHMPIDDigit::SizeAllY()+5); 
114   TLatex txt; txt.SetTextSize(0.1); txt.DrawLatex(-5,-5,Form("%i",iCh));
115   
116   for(Int_t iPc=AliHMPIDDigit::kMinPc;iPc<=AliHMPIDDigit::kMaxPc;iPc++){
117     TBox *pBox=new TBox(AliHMPIDDigit::fMinPcX[iPc],AliHMPIDDigit::fMinPcY[iPc],
118                         AliHMPIDDigit::fMaxPcX[iPc],AliHMPIDDigit::fMaxPcY[iPc]);
119     
120     pBox->SetFillStyle(0);
121     pBox->Draw();
122     
123   }//PC loop      
124 }//DrawCh()
125 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
126 void DrawEvt(TClonesArray *pHitLst,TObjArray *pDigLst,TObjArray *pCluLst,AliESD *pEsd)
127 {//draws all the objects of current event in given canvas
128
129   AliHMPIDRecon rec;  
130   TPolyMarker *pTxC[7];  TPolyMarker *pRin[7]; //intesections and rings
131   for(Int_t ch=0;ch<7;ch++){
132     pTxC[ch]=new TPolyMarker; pTxC[ch]->SetMarkerStyle(2); pTxC[ch]->SetMarkerColor(kRed); pTxC[ch]->SetMarkerSize(3);
133     pRin[ch]=new TPolyMarker; pRin[ch]->SetMarkerStyle(6); pRin[ch]->SetMarkerColor(kMagenta);
134   }
135   
136   
137   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop to collect cerenkov rings and intersection points
138     AliESDtrack *pTrk=pEsd->GetTrack(iTrk);
139     Int_t ch=pTrk->GetHMPIDcluIdx();
140     if(ch<0) continue; //this track does not hit HMPID
141     ch/=1000000; 
142     Float_t th,ph,xPc,yPc; pTrk->GetHMPIDtrk(xPc,yPc,th,ph);  //get info on current track
143     pTxC[ch]->SetNextPoint(xPc,yPc);                          //add new intersection point
144     
145     Float_t ckov=pTrk->GetHMPIDsignal();  Float_t err=TMath::Sqrt(pTrk->GetHMPIDchi2());
146     
147     if(ckov>0){
148       rec.SetTrack(xPc,yPc,th,ph);
149      TVector2 pos;  for(int j=0;j<100;j++){rec.TracePhot(ckov,j*0.0628,pos); pRin[ch]->SetNextPoint(pos.X(),pos.Y());}      
150     }
151   }//tracks loop
152       
153   for(Int_t iCh=0;iCh<7;iCh++){//chambers loop    
154     switch(iCh){
155       case 6: pAll->cd(1); break; case 5: pAll->cd(2); break;
156       case 4: pAll->cd(4); break; case 3: pAll->cd(5); break; case 2: pAll->cd(6); break;
157                                   case 1: pAll->cd(8); break; case 0: pAll->cd(9); break;
158     }
159     gPad->SetEditable(kTRUE); gPad->Clear(); 
160     DrawCh(iCh);
161     for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){
162       AliHMPIDHit *pHit=(AliHMPIDHit*)pHitLst->At(iHit);
163       if(pHit->Ch()==iCh)        pHit->Draw();  //draw hits
164     }
165     ((TClonesArray*)pDigLst->At(iCh))->Draw();  //draw digits
166     ((TClonesArray*)pCluLst->At(iCh))->Draw();  //draw clusters
167                             pTxC[iCh]->Draw();  //draw intersections
168                             pRin[iCh]->Draw();  //draw rings
169     gPad->SetEditable(kFALSE);
170   }//chambers loop
171 }//DrawEvt()
172 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
173 void DoZoom(Int_t evt, Int_t px, Int_t py, TObject *obj)
174 {
175   if(evt!=5 && evt!=6) return; //5- zoom in 6-zoom out
176   const Int_t minZoom=64;
177   const Int_t maxZoom=2;
178   static Int_t zoom=minZoom; //zoom level
179   if(evt==5&&zoom==maxZoom) return; 
180   if(evt==6&&zoom==minZoom) return; 
181   
182   if(!obj->IsA()->InheritsFrom("TPad")) return;  //current object is not pad
183   TPad *pPad=(TPad*)obj;
184   if(pPad->GetNumber()==3 || pPad->GetNumber()==7) return; //current pad is wrong
185
186  // Printf("evt=%i (%i,%i) %s",evt,px,py,obj->GetName());
187     
188   Float_t x=pPad->AbsPixeltoX(px); Float_t y=pPad->AbsPixeltoY(py); 
189  
190   if(evt==5){ zoom=zoom/2;     pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom in
191   else      { zoom=zoom*2;     pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom out 
192   if(zoom==minZoom) pPad->Range(-10,-10,AliHMPIDDigit::SizeAllX()+5,AliHMPIDDigit::SizeAllY()+5);
193   ((TCanvas *)gTQSender)->SetTitle(Form("zoom x%i",minZoom/zoom));
194   pPad->Modified();
195   pPad->Update();                                              
196 }
197 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++