]> git.uio.no Git - u/mrichter/AliRoot.git/blame - RICH/AliRICHDisplFast.cxx
TFlukaCerenkov added.
[u/mrichter/AliRoot.git] / RICH / AliRICHDisplFast.cxx
CommitLineData
6b17b320 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 **************************************************************************/
6b17b320 15#include "AliRICHDisplFast.h"
53fd478b 16#include "AliRICH.h"
17#include "AliRICHChamber.h"
6b17b320 18#include "AliRICHParam.h"
19#include <AliLoader.h>
6b17b320 20#include <TCanvas.h>
af3d25a6 21#include <TPolyLine.h>
6b17b320 22#include <TParticle.h>
23#include <TStyle.h>
6b17b320 24#include <TH2.h>
25#include <TMath.h>
d4c94996 26#include <TLatex.h>
6b17b320 27
aed240d4 28ClassImp(AliRICHDisplFast)
6b17b320 29
aed240d4 30//__________________________________________________________________________________________________
31void AliRICHDisplFast::Exec()
6b17b320 32{
d4c94996 33 AliRICH *pRich = (AliRICH*)gAlice->GetDetector("RICH");
34 Bool_t isHits =!pRich->GetLoader()->LoadHits();
35 Bool_t isDigits =!pRich->GetLoader()->LoadDigits();
36 Bool_t isClusters=!pRich->GetLoader()->LoadRecPoints();
37
38 if(!isHits && !isDigits && !isClusters){Error("Exec","No hits digits and clusters. Nothing to display.");return;}
39
40 TCanvas *Display = new TCanvas("Display","RICH Display",0,0,600,600);
41 gStyle->SetPalette(1);
42
43 TH2F *pHitsH2=0,*pDigitsH2=0,*pClustersH2=0;
44
45 if(isHits) pHitsH2 = new TH2F("pHitsH2" , "Event Display",165,-AliRICHParam::PcSizeX()/2,AliRICHParam::PcSizeX()/2,
af3d25a6 46 144,-AliRICHParam::PcSizeY()/2,AliRICHParam::PcSizeY()/2);
d4c94996 47 if(pHitsH2) pHitsH2->SetStats(kFALSE);
48
49 if(isDigits) pDigitsH2 = new TH2F("pDigitsH2" ,"Event Display",165,-AliRICHParam::PcSizeX()/2,AliRICHParam::PcSizeX()/2,
af3d25a6 50 144,-AliRICHParam::PcSizeY()/2,AliRICHParam::PcSizeY()/2);
d4c94996 51 if(isClusters) pClustersH2 = new TH2F("pClustersH2","Event Display",165,-AliRICHParam::PcSizeX()/2,AliRICHParam::PcSizeX()/2,
af3d25a6 52 144,-AliRICHParam::PcSizeY()/2,AliRICHParam::PcSizeY()/2);
aed240d4 53
6b17b320 54
c94cc386 55
d4c94996 56
aed240d4 57 for(Int_t iEventN=0;iEventN<gAlice->GetEventsPerRun();iEventN++){//events Loop
58 pRich->GetLoader()->GetRunLoader()->GetEvent(iEventN);
6b17b320 59
c94cc386 60
aed240d4 61 Int_t nPrimaries = (Int_t)pRich->GetLoader()->TreeH()->GetEntries();
62 TObjArray * Hits = new TObjArray[nPrimaries];
6b17b320 63
aed240d4 64 for(Int_t i=0;i<nPrimaries;i++) {
65 pRich->GetLoader()->TreeH()->GetEntry(i);
66 Int_t nHits = pRich->Hits()->GetEntries();
67 for(Int_t k=0;k<nHits;k++) Hits[i].Add(pRich->Hits()->At(k));
68
6b17b320 69 }
d4c94996 70//display all the staff on chamber by chamber basis
71 for(Int_t iChamber=1;iChamber<=7;iChamber++){//chambers loop
72 if(isHits) pHitsH2 ->Reset();
73 if(isDigits) pDigitsH2 ->Reset();
74 if(isClusters) pClustersH2->Reset();
75//deals with hits
aed240d4 76 for(Int_t i=0;i<nPrimaries;i++){//prims loop
6b17b320 77 pRich->GetLoader()->TreeH()->GetEntry(i);
78 Int_t nHits = pRich->Hits()->GetEntries();
aed240d4 79 for(Int_t j=0;j<nHits;j++){//hits loop
6b17b320 80 AliRICHhit *pHit = (AliRICHhit*)Hits[i].At(j);
aed240d4 81 if(pHit->C()==iChamber){
d4c94996 82 TVector3 hitGlobX3= pHit->OutX3();
83 TVector2 hitLocX2 = pRich->C(iChamber)->Glob2Loc(hitGlobX3);
84 pHitsH2->Fill(hitLocX2.X(),hitLocX2.Y(),200);
aed240d4 85 }//if
86 }//hits loop
87 }//prims loop
d4c94996 88 pHitsH2->SetTitle(Form("event %i chamber %2i",iEventN,iChamber));
aed240d4 89 pHitsH2->SetMarkerColor(kRed); pHitsH2->SetMarkerStyle(29); pHitsH2->SetMarkerSize(0.4);
90 pHitsH2->Draw();
af3d25a6 91 DrawSectors();
d4c94996 92 TLatex l; l.SetNDC(); l.SetTextSize(0.02);
93 if(!isHits) {l.SetTextColor(kRed) ;l.DrawLatex(0.1,0.01,"No Hits" );}
94 if(!isDigits) {l.SetTextColor(kGreen);l.DrawLatex(0.4,0.01,"No DIGITS" );}
95 if(!isClusters) {l.SetTextColor(kBlue) ;l.DrawLatex(0.8,0.01,"No CLUSTERS");}
6b17b320 96 Display->Update();
97 Display->Modified();
d4c94996 98 getchar();
99//deals with digits
100 if(isDigits){
101 pRich->GetLoader()->TreeD()->GetEntry(0);
102 for(Int_t j=0;j<pRich->Digits(iChamber)->GetEntries();j++){//digits loop
103 AliRICHdigit *pDig = (AliRICHdigit*)pRich->Digits(iChamber)->At(j);
104 TVector2 x2=AliRICHParam::Pad2Loc(pDig->X(),pDig->Y());
105 pDigitsH2->Fill(x2.X(),x2.Y(),100);
106 }//digits loop
107 pDigitsH2->SetMarkerColor(kGreen); pDigitsH2->SetMarkerStyle(29); pDigitsH2->SetMarkerSize(0.4);
108 pDigitsH2->Draw("same");
109 Display->Update();
110 Display->Modified();
111 getchar();
112 }//if(isDigits)
113//deals with clusters
114 if(isClusters){
115 pRich->GetLoader()->TreeR()->GetEntry(0);
116 for(Int_t j=0;j<pRich->Clusters(iChamber)->GetEntries();j++){//clusters loop
117 AliRICHcluster *pClus = (AliRICHcluster*)pRich->Clusters(iChamber)->At(j);
118 pClustersH2->Fill(pClus->X(),pClus->Y(),50);
119 }//clusters loop
120 pClustersH2->SetMarkerColor(kBlue); pClustersH2->SetMarkerStyle(29); pClustersH2->SetMarkerSize(0.4);
121 pClustersH2->Draw("same");
122 Display->Update();
123 Display->Modified();
124 getchar();
125 }//if(isClusters)
126 }//chambers loop
aed240d4 127 delete [] Hits;
d4c94996 128 }//events Loop
aed240d4 129 pRich->GetLoader()->UnloadHits();
d4c94996 130 if(isDigits) pRich->GetLoader()->UnloadDigits();
131 if(isClusters) pRich->GetLoader()->UnloadRecPoints();
aed240d4 132}//Exec()
133//__________________________________________________________________________________________________
af3d25a6 134void AliRICHDisplFast::DrawSectors()
135{
136 Double_t x1[5] = {-AliRICHParam::PcSizeX()/2,
137 -AliRICHParam::SectorSizeX()/2-AliRICHParam::DeadZone(),
138 -AliRICHParam::SectorSizeX()/2-AliRICHParam::DeadZone(),
139 -AliRICHParam::PcSizeX()/2,
140 -AliRICHParam::PcSizeX()/2};
141 Double_t y1[5] = {AliRICHParam::DeadZone()/2,
142 AliRICHParam::DeadZone()/2,
143 AliRICHParam::PcSizeY()/2,
144 AliRICHParam::PcSizeY()/2,
145 AliRICHParam::DeadZone()/2};
146 Double_t x2[5] = {-AliRICHParam::SectorSizeX()/2,
147 AliRICHParam::SectorSizeX()/2,
148 AliRICHParam::SectorSizeX()/2,
149 -AliRICHParam::SectorSizeX()/2,
150 -AliRICHParam::SectorSizeX()/2};
151 Double_t y2[5] = {AliRICHParam::DeadZone()/2,
152 AliRICHParam::DeadZone()/2,
153 AliRICHParam::PcSizeY()/2,
154 AliRICHParam::PcSizeY()/2,
155 AliRICHParam::DeadZone()/2};
156 Double_t x3[5] = { AliRICHParam::SectorSizeX()/2+AliRICHParam::DeadZone(),
157 AliRICHParam::PcSizeX()/2,
158 AliRICHParam::PcSizeX()/2,
159 AliRICHParam::SectorSizeX()/2+AliRICHParam::DeadZone(),
160 AliRICHParam::SectorSizeX()/2+AliRICHParam::DeadZone()};
161 Double_t y3[5] = {AliRICHParam::DeadZone()/2,
162 AliRICHParam::DeadZone()/2,
163 AliRICHParam::PcSizeY()/2,
164 AliRICHParam::PcSizeY()/2,
165 AliRICHParam::DeadZone()/2};
166 Double_t x4[5] = {-AliRICHParam::PcSizeX()/2,
167 -AliRICHParam::SectorSizeX()/2-AliRICHParam::DeadZone(),
168 -AliRICHParam::SectorSizeX()/2-AliRICHParam::DeadZone(),
169 -AliRICHParam::PcSizeX()/2,
170 -AliRICHParam::PcSizeX()/2};
171 Double_t y4[5] = {-AliRICHParam::PcSizeY()/2,
172 -AliRICHParam::PcSizeY()/2,
173 -AliRICHParam::DeadZone()/2,
174 -AliRICHParam::DeadZone()/2,
175 -AliRICHParam::PcSizeY()/2};
176 Double_t x5[5] = {-AliRICHParam::SectorSizeX()/2,
177 AliRICHParam::SectorSizeX()/2,
178 AliRICHParam::SectorSizeX()/2,
179 -AliRICHParam::SectorSizeX()/2,
180 -AliRICHParam::SectorSizeX()/2};
181 Double_t y5[5] = {-AliRICHParam::PcSizeY()/2,
182 -AliRICHParam::PcSizeY()/2,
183 -AliRICHParam::DeadZone()/2,
184 -AliRICHParam::DeadZone()/2,
185 -AliRICHParam::PcSizeY()/2};
186 Double_t x6[5] = { AliRICHParam::SectorSizeX()/2+AliRICHParam::DeadZone(),
187 AliRICHParam::PcSizeX()/2,
188 AliRICHParam::PcSizeX()/2,
189 AliRICHParam::SectorSizeX()/2+AliRICHParam::DeadZone(),
190 AliRICHParam::SectorSizeX()/2+AliRICHParam::DeadZone()};
191 Double_t y6[5] = {-AliRICHParam::PcSizeY()/2,
192 -AliRICHParam::PcSizeY()/2,
193 -AliRICHParam::DeadZone()/2,
194 -AliRICHParam::DeadZone()/2,
195 -AliRICHParam::PcSizeY()/2};
d4c94996 196 TPolyLine *sector1 = new TPolyLine(5,x1,y1); sector1->SetLineColor(21); sector1->Draw();
197 TPolyLine *sector2 = new TPolyLine(5,x2,y2); sector2->SetLineColor(21); sector2->Draw();
198 TPolyLine *sector3 = new TPolyLine(5,x3,y3); sector3->SetLineColor(21); sector3->Draw();
199 TPolyLine *sector4 = new TPolyLine(5,x4,y4); sector4->SetLineColor(21); sector4->Draw();
200 TPolyLine *sector5 = new TPolyLine(5,x5,y5); sector5->SetLineColor(21); sector5->Draw();
201 TPolyLine *sector6 = new TPolyLine(5,x6,y6); sector6->SetLineColor(21); sector6->Draw();
202}//DrawSectors()