]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHDisplFast.cxx
No return in void function (Alpha)
[u/mrichter/AliRoot.git] / RICH / AliRICHDisplFast.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 "AliRICHDisplFast.h"
16 #include "AliRICH.h"
17 #include "AliRICHChamber.h"
18 #include "AliRICHParam.h"
19 #include <AliLoader.h>
20 #include <TCanvas.h>
21 #include <TPolyLine.h>
22 #include <TParticle.h>
23 #include <TStyle.h>
24 #include <TMath.h>
25 #include <TLatex.h>
26
27 ClassImp(AliRICHDisplFast)
28 //__________________________________________________________________________________________________
29 void AliRICHDisplFast::ShowEvent(Int_t iEvtNmin,Int_t iEvtNmax)
30 {
31   TH2F *pDigitsH2[8];
32   char titobj[11],titdisp[20];
33
34   AliRICH *pRich = (AliRICH*)gAlice->GetDetector("RICH");
35   Bool_t isDigits  =!pRich->GetLoader()->LoadDigits();
36   if(!isDigits){Error("ShoEvent","No digits. Nothing to display.");return;}
37   
38   TCanvas *canvas = new TCanvas("RICHDisplay","RICH Display",0,0,1226,900);
39 //  TCanvas *canvas = (TCanvas*)gROOT->FindObjectAny("RICHDisplay");
40   
41   gStyle->SetPalette(1);
42
43   canvas->Divide(3,3);
44   
45   for(Int_t iChamber=1;iChamber<=7;iChamber++) {
46     sprintf(titobj,"pDigitsH2_%i",iChamber);
47     sprintf(titdisp,"Chamber  %i",iChamber);
48     pDigitsH2[iChamber] = new TH2F(titobj,titdisp,165,0,AliRICHParam::PcSizeX(),144,0,AliRICHParam::PcSizeY());
49     pDigitsH2[iChamber]->SetMarkerColor(kGreen); 
50     pDigitsH2[iChamber]->SetMarkerStyle(29); 
51     pDigitsH2[iChamber]->SetMarkerSize(0.4);
52     pDigitsH2[iChamber]->SetStats(kFALSE);
53   }
54   
55   if(iEvtNmax>gAlice->GetEventsPerRun()) iEvtNmax=gAlice->GetEventsPerRun();
56
57   TText *t = new TText();
58   t->SetTextSize(0.10);
59
60   TText *tit = new TText(0.1,0.6,"RICH Display");
61   tit->SetTextSize(0.10);
62
63   for(Int_t iEventN=iEvtNmin;iEventN<=iEvtNmax;iEventN++) {
64         
65     canvas->cd(1);
66     sprintf(titdisp,"Event Number %i",iEventN);
67     t->SetText(0.2,0.4,titdisp);
68     t->Draw();
69     tit->Draw();
70     pRich->GetLoader()->GetRunLoader()->GetEvent(iEventN);
71     pRich->GetLoader()->TreeD()->GetEntry(0);
72     for(Int_t iChamber=1;iChamber<=7;iChamber++) {
73       pDigitsH2[iChamber]->Reset();    
74       for(Int_t j=0;j<pRich->Digits(iChamber)->GetEntries();j++) {//digits loop
75         AliRICHdigit *pDig = (AliRICHdigit*)pRich->Digits(iChamber)->At(j);
76         TVector2 x2=AliRICHParam::Pad2Loc(pDig->Pad());
77         pDigitsH2[iChamber]->Fill(x2.X(),x2.Y());
78       }
79       if(iChamber==1) canvas->cd(7);
80       if(iChamber==2) canvas->cd(8);
81       if(iChamber==3) canvas->cd(4);
82       if(iChamber==4) canvas->cd(5);
83       if(iChamber==5) canvas->cd(6);
84       if(iChamber==6) canvas->cd(2);
85       if(iChamber==7) canvas->cd(3);
86       pDigitsH2[iChamber]->Draw();
87     }
88     canvas->Update();
89     canvas->Modified();
90     
91     if(iEvtNmin<iEvtNmax) gPad->WaitPrimitive();
92 //    canvas->cd(3);
93 //    gPad->Clear();
94   }
95 }
96 //__________________________________________________________________________________________________
97 void AliRICHDisplFast::Exec(Option_t *)
98 {
99   AliRICH *pRich = (AliRICH*)gAlice->GetDetector("RICH");
100   Bool_t isHits    =!pRich->GetLoader()->LoadHits();
101   Bool_t isDigits  =!pRich->GetLoader()->LoadDigits();
102   Bool_t isClusters=!pRich->GetLoader()->LoadRecPoints();
103   
104   if(!isHits && !isDigits && !isClusters){Error("Exec","No hits digits and clusters. Nothing to display.");return;}
105   
106   TCanvas *Display = new TCanvas("Display","RICH Display",0,0,600,600);
107   gStyle->SetPalette(1);
108   
109   TH2F *pHitsH2=0,*pDigitsH2=0,*pClustersH2=0;
110   
111   if(isHits)     pHitsH2     = new TH2F("pHitsH2"  ,  "Event Display;x,cm;y,cm",165,0,AliRICHParam::PcSizeX(),
112                                                                       144,0,AliRICHParam::PcSizeY());
113   if(pHitsH2)    pHitsH2->SetStats(kFALSE);
114   
115   if(isDigits)   pDigitsH2   = new TH2F("pDigitsH2"  ,"Event Display",165,0,AliRICHParam::PcSizeX(),
116                                                                       144,0,AliRICHParam::PcSizeY());
117   if(isClusters) pClustersH2 = new TH2F("pClustersH2","Event Display",165,0,AliRICHParam::PcSizeX(),
118                                                                       144,0,AliRICHParam::PcSizeY());
119   
120   for(Int_t iEventN=0;iEventN<gAlice->GetEventsPerRun();iEventN++){//events Loop
121     pRich->GetLoader()->GetRunLoader()->GetEvent(iEventN);
122
123   
124     Int_t nPrimaries = (Int_t)pRich->GetLoader()->TreeH()->GetEntries();
125     TObjArray * Hits = new TObjArray[nPrimaries];
126   
127     for(Int_t i=0;i<nPrimaries;i++) {
128       pRich->GetLoader()->TreeH()->GetEntry(i);
129       Int_t nHits = pRich->Hits()->GetEntries();
130       for(Int_t k=0;k<nHits;k++)         Hits[i].Add(pRich->Hits()->At(k));
131     
132     }
133 //display all the staff on chamber by chamber basis           
134     for(Int_t iChamber=1;iChamber<=7;iChamber++){//chambers loop       
135       if(isHits)     pHitsH2    ->Reset();     
136       if(isDigits)   pDigitsH2  ->Reset();     
137       if(isClusters) pClustersH2->Reset();
138 //deals with hits
139       for(Int_t i=0;i<nPrimaries;i++){//prims loop
140         pRich->GetLoader()->TreeH()->GetEntry(i);
141         Int_t nHits = pRich->Hits()->GetEntries();
142         for(Int_t j=0;j<nHits;j++){//hits loop
143           AliRICHhit *pHit = (AliRICHhit*)Hits[i].At(j);
144           if(pHit->C()==iChamber){
145             TVector3 hitGlobX3= pHit->OutX3();
146             TVector2 hitLocX2 = pRich->C(iChamber)->Mrs2Pc(hitGlobX3);
147             pHitsH2->Fill(hitLocX2.X(),hitLocX2.Y(),200);
148           }//if
149         }//hits loop         
150       }//prims loop
151       pHitsH2->SetTitle(Form("event %i chamber %2i",iEventN,iChamber));
152       pHitsH2->SetMarkerColor(kRed); pHitsH2->SetMarkerStyle(29); pHitsH2->SetMarkerSize(0.4);
153       pHitsH2->Draw();
154       DrawSectors();
155       TLatex l; l.SetNDC(); l.SetTextSize(0.02);
156       if(!isHits)     {l.SetTextColor(kRed)  ;l.DrawLatex(0.1,0.01,"No Hits"    );}
157       if(!isDigits)   {l.SetTextColor(kGreen);l.DrawLatex(0.4,0.01,"No DIGITS"  );}
158       if(!isClusters) {l.SetTextColor(kBlue) ;l.DrawLatex(0.8,0.01,"No CLUSTERS");}
159       Display->Update();
160       Display->Modified();
161       gPad->WaitPrimitive();
162 //deals with digits      
163       if(isDigits){
164         pRich->GetLoader()->TreeD()->GetEntry(0);
165         for(Int_t j=0;j<pRich->Digits(iChamber)->GetEntries();j++){//digits loop
166           AliRICHdigit *pDig = (AliRICHdigit*)pRich->Digits(iChamber)->At(j);
167           TVector2 x2=AliRICHParam::Pad2Loc(pDig->Pad());
168           pDigitsH2->Fill(x2.X(),x2.Y(),100);
169         }//digits loop
170         pDigitsH2->SetMarkerColor(kGreen); pDigitsH2->SetMarkerStyle(29); pDigitsH2->SetMarkerSize(0.4);
171         pDigitsH2->Draw("same");
172         Display->Update();
173         Display->Modified();       
174         gPad->WaitPrimitive();
175       }//if(isDigits)      
176 //deals with clusters      
177       if(isClusters){
178         pRich->GetLoader()->TreeR()->GetEntry(0);
179         for(Int_t j=0;j<pRich->Clusters(iChamber)->GetEntries();j++){//clusters loop
180           AliRICHcluster *pClus = (AliRICHcluster*)pRich->Clusters(iChamber)->At(j);
181           pClustersH2->Fill(pClus->X(),pClus->Y(),50);
182         }//clusters loop
183         pClustersH2->SetMarkerColor(kBlue); pClustersH2->SetMarkerStyle(29);  pClustersH2->SetMarkerSize(0.4);
184         pClustersH2->Draw("same");
185         Display->Update();
186         Display->Modified();
187         gPad->WaitPrimitive();
188       }//if(isClusters)
189     }//chambers loop
190     delete [] Hits;
191   }//events Loop
192   pRich->GetLoader()->UnloadHits();
193   if(isDigits)   pRich->GetLoader()->UnloadDigits();
194   if(isClusters) pRich->GetLoader()->UnloadRecPoints();
195 }//Exec()
196 //__________________________________________________________________________________________________
197 void AliRICHDisplFast::DrawSectors() 
198
199   AliRICHParam p;
200   Double_t xLeft[5]  = {0,0,p.SectorSizeX(),p.SectorSizeX(),0};
201   Double_t xRight[5] = {p.SectorSizeX()+p.DeadZone(),p.SectorSizeX()+p.DeadZone(),p.PcSizeX(),p.PcSizeX(),p.SectorSizeX()+p.DeadZone()};
202   
203   Double_t yDown[5]   = {0,p.SectorSizeY(),p.SectorSizeY(),0,0};
204   Double_t yCenter[5] = {  p.SectorSizeY()+p.DeadZone(),2*p.SectorSizeY()+p.DeadZone(),2*p.SectorSizeY()+p.DeadZone(),
205                            p.SectorSizeY()+p.DeadZone(),p.SectorSizeY()+p.DeadZone()};  
206   Double_t yUp[5]     = {2*p.SectorSizeY()+2*p.DeadZone(),p.PcSizeY(),p.PcSizeY(),2*p.SectorSizeY()+2*p.DeadZone(),2*p.SectorSizeY()+2*p.DeadZone()};
207   
208   TPolyLine *sec1 = new TPolyLine(5,xLeft ,yDown);    sec1->SetLineColor(21);  sec1->Draw();
209   TPolyLine *sec2 = new TPolyLine(5,xRight,yDown);    sec2->SetLineColor(21);  sec2->Draw();
210   TPolyLine *sec3 = new TPolyLine(5,xLeft ,yCenter);  sec3->SetLineColor(21);  sec3->Draw();
211   TPolyLine *sec4 = new TPolyLine(5,xRight,yCenter);  sec4->SetLineColor(21);  sec4->Draw();
212   TPolyLine *sec5 = new TPolyLine(5,xLeft, yUp);      sec5->SetLineColor(21);  sec5->Draw();
213   TPolyLine *sec6 = new TPolyLine(5,xRight,yUp);      sec6->SetLineColor(21);  sec6->Draw();
214 }//DrawSectors()