]>
Commit | Line | Data |
---|---|---|
c0d7adf8 | 1 | #if !defined(__CINT__) || defined(__MAKECINT__) |
2 | ||
e9f11028 | 3 | #include <TSystem.h> |
4 | #include <TFile.h> | |
5 | #include <TTree.h> | |
6 | #include <TButton.h> | |
7 | #include <TCanvas.h> | |
45b05be5 | 8 | #include <TString.h> |
e9f11028 | 9 | #include <TClonesArray.h> |
10 | #include <TParticle.h> | |
11 | #include <TRandom.h> | |
12 | #include <TLatex.h> | |
45b05be5 | 13 | #include <TPDGCode.h> |
e9f11028 | 14 | #include <TLegend.h> |
15 | #include <TPolyMarker.h> | |
16 | #include <TBox.h> | |
6e06db1d | 17 | #include <AliESDEvent.h> |
e9f11028 | 18 | #include <AliCDBManager.h> |
19 | #include <AliCDBEntry.h> | |
20 | #include "AliHMPIDHit.h" | |
21 | #include "AliHMPIDv2.h" | |
22 | #include "AliHMPIDReconstructor.h" | |
23 | #include "AliHMPIDRecon.h" | |
24 | #include "AliHMPIDParam.h" | |
25 | #include "AliHMPIDCluster.h" | |
e4290ede | 26 | #include "AliTracker.h" |
27 | ||
e9f11028 | 28 | #include <TChain.h> |
126b3e0e | 29 | |
c0d7adf8 | 30 | #endif |
31 | ||
45b05be5 | 32 | TCanvas *fCanvas=0; Int_t fType=3; Int_t fEvt=-1; Int_t fNevt=0; |
e9f11028 | 33 | TFile *fHitFile; TTree *fHitTree; TClonesArray *fHitLst; TPolyMarker *fRenMip[7]; TPolyMarker *fRenCko[7]; TPolyMarker *fRenFee[7]; |
34 | TClonesArray *fSdiLst; | |
35 | TFile *fDigFile; TTree *fDigTree; TObjArray *fDigLst; TPolyMarker *fRenDig[7]; | |
36 | TFile *fCluFile; TTree *fCluTree; TObjArray *fCluLst; TPolyMarker *fRenClu[7]; | |
45b05be5 | 37 | TFile *fEsdFile; TTree *fEsdTree; AliESDEvent *fEsd; TPolyMarker *fRenTxC[7]; TPolyMarker *fRenRin[7]; |
e9f11028 | 38 | TFile *fCosFile; TTree *fCosTree; |
45b05be5 | 39 | |
40 | TButton *fHitMipBok,*fHitCkoBok,*fHitFeeBok,*fDigBok,*fCluBok,*fEsdBok; | |
41 | ||
42 | TString fStHitMip = "ON"; | |
43 | TString fStHitCko = "ON"; | |
44 | TString fStHitFee = "ON"; | |
45 | TString fStDig = "ON"; | |
46 | TString fStClu = "ON"; | |
47 | TString fStEsd = "ON"; | |
48 | ||
49 | Int_t fTimeArrow = 1; | |
50 | ||
e4290ede | 51 | AliRunLoader *gAL=0; |
8282f9d0 | 52 | |
e9f11028 | 53 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
54 | void CreateContainers() | |
55 | {//to create all containers | |
56 | fHitLst=new TClonesArray("AliHMPIDHit"); | |
57 | fSdiLst=new TClonesArray("AliHMPIDDigit"); | |
58 | fDigLst=new TObjArray(7); for(Int_t i=0;i<7;i++) fDigLst->AddAt(new TClonesArray("AliHMPIDDigit"),i); fDigLst->SetOwner(kTRUE); | |
59 | fCluLst=new TObjArray(7); for(Int_t i=0;i<7;i++) fCluLst->AddAt(new TClonesArray("AliHMPIDCluster"),i); fCluLst->SetOwner(kTRUE); | |
6e06db1d | 60 | fEsd =new AliESDEvent; |
126b3e0e | 61 | } |
126b3e0e | 62 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
e9f11028 | 63 | void CreateRenders() |
64 | { | |
65 | for(Int_t ch=0;ch<7;ch++){ | |
66 | fRenMip[ch]=new TPolyMarker; fRenMip[ch]->SetMarkerStyle(kOpenTriangleUp); fRenMip[ch]->SetMarkerColor(kRed); | |
67 | fRenCko[ch]=new TPolyMarker; fRenCko[ch]->SetMarkerStyle(kOpenCircle); fRenCko[ch]->SetMarkerColor(kRed); | |
68 | fRenFee[ch]=new TPolyMarker; fRenFee[ch]->SetMarkerStyle(kOpenDiamond); fRenFee[ch]->SetMarkerColor(kRed); | |
45b05be5 | 69 | fRenDig[ch]=new TPolyMarker; fRenDig[ch]->SetMarkerStyle(kOpenSquare); fRenDig[ch]->SetMarkerColor(kGreen); |
e9f11028 | 70 | fRenClu[ch]=new TPolyMarker; fRenClu[ch]->SetMarkerStyle(kStar); fRenClu[ch]->SetMarkerColor(kBlue); |
71 | fRenTxC[ch]=new TPolyMarker; fRenTxC[ch]->SetMarkerStyle(kPlus); fRenTxC[ch]->SetMarkerColor(kRed); fRenTxC[ch]->SetMarkerSize(3); | |
72 | fRenRin[ch]=new TPolyMarker; fRenRin[ch]->SetMarkerStyle(kFullDotSmall); fRenRin[ch]->SetMarkerColor(kMagenta); | |
73 | } | |
74 | }//CreateRenders() | |
75 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
76 | void ClearRenders() | |
77 | { | |
78 | for(Int_t ch=0;ch<7;ch++){ | |
79 | fRenTxC[ch]->SetPolyMarker(0); | |
80 | fRenRin[ch]->SetPolyMarker(0); | |
81 | fRenMip[ch]->SetPolyMarker(0); | |
82 | fRenCko[ch]->SetPolyMarker(0); | |
83 | fRenFee[ch]->SetPolyMarker(0); | |
84 | fRenDig[ch]->SetPolyMarker(0); | |
85 | fRenClu[ch]->SetPolyMarker(0); | |
86 | } | |
87 | }//ClearRenders() | |
88 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
89 | void PrintHits() | |
90 | { | |
91 | //Prints a list of HMPID hits for a given event. Default is event number 0. | |
92 | if(!fHitTree) return; | |
93 | Printf("List of HMPID hits for event %i",fEvt); | |
94 | Int_t iTot=0; | |
95 | for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){//entries loop | |
96 | fHitTree->GetEntry(iEnt); | |
97 | fHitLst->Print(); | |
98 | iTot+=fHitLst->GetEntries(); | |
99 | } | |
100 | Printf("totally %i hits for event %i",iTot,fEvt); | |
101 | }//PrintHits(); | |
102 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
103 | void PrintSdis() | |
104 | {//prints a list of HMPID sdigits for a given event | |
105 | Printf("List of HMPID sdigits for event %i",fEvt); | |
106 | fSdiLst->Print(); | |
107 | Printf("totally %i sdigits for event %i",fSdiLst->GetEntries(),fEvt); | |
108 | }//PrintSdis() | |
109 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
110 | void PrintDigs() | |
111 | {//prints a list of HMPID digits | |
112 | Printf("List of HMPID digits for event %i",fEvt); | |
113 | fDigLst->Print(); | |
114 | Int_t iTot=0; for(Int_t i=0;i<7;i++) {iTot+=((TClonesArray*)fDigLst->At(i))->GetEntries();} | |
115 | Printf("totally %i digits for event %i",iTot,fEvt); | |
8282f9d0 | 116 | } |
117 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
e9f11028 | 118 | void PrintClus() |
119 | {//prints a list of HMPID clusters for a given event | |
120 | Printf("List of HMPID clusters for event %i",fEvt); | |
8282f9d0 | 121 | |
e9f11028 | 122 | fCluLst->Print(); |
123 | Int_t iTot=0; for(Int_t iCh=0;iCh<7;iCh++) iTot+=((TClonesArray*)fCluLst->At(iCh))->GetEntries(); | |
124 | ||
125 | Printf("totally %i clusters for event %i",iTot,fEvt); | |
8282f9d0 | 126 | } |
127 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
e9f11028 | 128 | void PrintEsd() |
129 | {//prints a list of HMPID Esd for a given event | |
130 | Printf("List of HMPID ESD summary for event %i",fEvt); | |
131 | for(Int_t iTrk=0;iTrk<fEsd->GetNumberOfTracks();iTrk++){ | |
132 | AliESDtrack *pTrk = fEsd->GetTrack(iTrk); | |
133 | Float_t x,y;Int_t q,nacc; pTrk->GetHMPIDmip(x,y,q,nacc); | |
134 | Printf("Track %02i with p %7.2f with ThetaCer %5.3f with %i photons",iTrk,pTrk->GetP(),pTrk->GetHMPIDsignal(),nacc); | |
135 | } | |
136 | }//PrintEsd() | |
137 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
138 | void DrawChamber(Int_t iCh) | |
139 | {//used by Draw() to Draw() chamber structure | |
45b05be5 | 140 | gPad->Range(-10,-10,AliHMPIDParam::SizeAllX()+5,AliHMPIDParam::SizeAllY()+5); |
e9f11028 | 141 | if(iCh>=0){TLatex txt; txt.SetTextSize(0.1); txt.DrawLatex(-5,-5,Form("%i",iCh));} |
142 | ||
ae5a42aa | 143 | for(Int_t iPc=AliHMPIDParam::kMinPc;iPc<=AliHMPIDParam::kMaxPc;iPc++){ |
144 | TBox *pBox=new TBox(AliHMPIDParam::MinPcX(iPc),AliHMPIDParam::MinPcY(iPc), | |
145 | AliHMPIDParam::MaxPcX(iPc),AliHMPIDParam::MaxPcY(iPc)); | |
e9f11028 | 146 | pBox->SetFillStyle(0); pBox->Draw(); |
147 | }//PC loop | |
148 | }//DrawChamber() | |
149 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
150 | void DrawLegend() | |
151 | {//used by Draw() to draw legend | |
152 | Int_t nTxC=0,nMip=0,nCko=0,nFee=0,nDig=0,nClu=0; | |
153 | for(Int_t ch=0;ch<7;ch++){ | |
45b05be5 | 154 | nTxC+=fRenTxC[ch]->GetLastPoint()+1; |
155 | nMip+=fRenMip[ch]->GetLastPoint()+1; | |
156 | nCko+=fRenCko[ch]->GetLastPoint()+1; | |
157 | nFee+=fRenFee[ch]->GetLastPoint()+1; | |
158 | nDig+=fRenDig[ch]->GetLastPoint()+1; | |
159 | nClu+=fRenClu[ch]->GetLastPoint()+1; | |
e9f11028 | 160 | } |
161 | TLegend *pLeg=new TLegend(0.2,0.2,0.8,0.8); | |
162 | pLeg->SetHeader(Form("Event %i Total %i",fEvt,fNevt)); | |
163 | pLeg->AddEntry(fRenTxC[0],Form("TRKxPC %i" ,nTxC),"p"); | |
164 | pLeg->AddEntry(fRenMip[0],Form("Mip hits %i" ,nMip),"p"); | |
165 | pLeg->AddEntry(fRenCko[0],Form("Ckov hits %i" ,nCko),"p"); | |
166 | pLeg->AddEntry(fRenFee[0],Form("Feed hits %i" ,nFee),"p"); | |
167 | pLeg->AddEntry(fRenDig[0],Form("Digs %i" ,nDig),"p"); | |
168 | pLeg->AddEntry(fRenClu[0],Form("Clus %i" ,nClu),"p"); | |
169 | pLeg->Draw(); | |
170 | }//DrawLegend() | |
171 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
172 | void Draw() | |
173 | {//draws all the objects of current event in given canvas | |
45b05be5 | 174 | Int_t iPadN[7]={9,8,6,5,4,2,1}; |
e9f11028 | 175 | for(Int_t iCh=0;iCh<7;iCh++){//chambers loop |
45b05be5 | 176 | fCanvas->cd(iPadN[iCh]); |
e9f11028 | 177 | gPad->SetEditable(kTRUE); gPad->Clear(); DrawChamber(iCh); |
45b05be5 | 178 | if(fStEsd =="ON") fRenTxC[iCh]->Draw(); |
179 | if(fStHitMip=="ON") fRenMip[iCh]->Draw(); | |
180 | if(fStHitFee=="ON") fRenFee[iCh]->Draw(); | |
181 | if(fStHitCko=="ON") fRenCko[iCh]->Draw(); | |
182 | if(fStEsd =="ON") fRenRin[iCh]->Draw(); | |
183 | if(fStDig =="ON") fRenDig[iCh]->Draw(); | |
184 | if(fStClu =="ON") fRenClu[iCh]->Draw(); | |
e9f11028 | 185 | gPad->SetEditable(kFALSE); |
186 | }//chambers loop | |
e9f11028 | 187 | fCanvas->cd(3); gPad->Clear(); DrawLegend(); |
45b05be5 | 188 | |
189 | fCanvas->cd(7); | |
190 | ||
191 | if(fHitFile){fHitMipBok->SetTitle(Form("Mips %s",fStHitMip.Data())); fHitMipBok->Modified();} | |
192 | if(fHitFile){fHitCkoBok->SetTitle(Form("Ckov %s",fStHitCko.Data())); fHitCkoBok->Modified();} | |
193 | if(fHitFile){fHitFeeBok->SetTitle(Form("Fdbk %s",fStHitFee.Data())); fHitFeeBok->Modified();} | |
194 | if(fDigFile){fDigBok->SetTitle(fStDig); fDigBok->Modified();} | |
195 | if(fCluFile){fCluBok->SetTitle(fStClu); fCluBok->Modified();} | |
196 | if(fEsdFile){fEsdBok->SetTitle(fStEsd); fEsdBok->Modified();} | |
197 | ||
e9f11028 | 198 | }//Draw() |
199 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
200 | void RenderHit(TClonesArray *pHitLst) | |
201 | {//used by ReadEvent() and SimulateEvent() to render hits to polymarker structures, one per chamber | |
202 | for(Int_t iHit=0;iHit<pHitLst->GetEntries();iHit++){ //hits loop | |
203 | AliHMPIDHit *pHit = (AliHMPIDHit*)pHitLst->At(iHit); Int_t ch=pHit->Ch(); Float_t x=pHit->LorsX(); Float_t y=pHit->LorsY(); //get current hit | |
204 | switch(pHit->Pid()){ | |
205 | case 50000050: fRenCko[ch]->SetNextPoint(x,y);break; | |
206 | case 50000051: fRenFee[ch]->SetNextPoint(x,y);break; | |
207 | default: fRenMip[ch]->SetNextPoint(x,y);break; | |
208 | }//switch hit PID | |
209 | }//hits loop for this entry | |
210 | }//RenderHits() | |
211 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
e4290ede | 212 | void RenderDig(TObjArray *pDigLst) |
213 | {//used by ReadEvent() and SimulateEvent() to render digs to Tbox structures, one per chamber | |
214 | for(Int_t iCh=0;iCh<=6;iCh++){ //chambers loop | |
215 | TClonesArray *pDigCham=(TClonesArray*)pDigLst->At(iCh); //get digs list for this chamber | |
216 | for(Int_t iDig=0;iDig<pDigCham->GetEntries();iDig++){ //digs loop | |
217 | AliHMPIDDigit *pDig = (AliHMPIDDigit*)pDigCham->At(iDig); Float_t x=pDig->LorsX(); Float_t y=pDig->LorsY(); //get current hit | |
218 | fRenDig[iCh]->SetNextPoint(x,y); | |
45b05be5 | 219 | }//Digit loop |
e4290ede | 220 | }//hits loop for this entry |
221 | }//RenderHits() | |
222 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
223 | void RenderClu(TObjArray *pClus) | |
e9f11028 | 224 | {//used by ReadEvent() and SimulateEvent() to render clusters to polymarker structures, one per chamber |
225 | for(Int_t iCh=0;iCh<=6;iCh++){ //chambers loop | |
226 | TClonesArray *pClusCham=(TClonesArray*)pClus->At(iCh); //get clusters list for this chamber | |
227 | for(Int_t iClu=0;iClu<pClusCham->GetEntries();iClu++){ //clusters loop | |
228 | AliHMPIDCluster *pClu = (AliHMPIDCluster*)pClusCham->At(iClu); //get current cluster | |
e4290ede | 229 | fRenClu[iCh]->SetNextPoint(pClu->X(),pClu->Y()); |
45b05be5 | 230 | // Printf("RenderClu: ch %i x %f y %f",iCh,pClu->X(),pClu->Y()); |
231 | }//cluster loop | |
232 | }//chamber loop for this entry | |
e9f11028 | 233 | }//RenderClus() |
234 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
6e06db1d | 235 | void RenderEsd(AliESDEvent *pEsd) |
e9f11028 | 236 | {//used by ReadEvent() or SimulateEvent() to render ESD to polymarker structures for rings and intersections one per chamber |
237 | AliHMPIDRecon rec; | |
238 | for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop to collect cerenkov rings and intersection points | |
239 | AliESDtrack *pTrk=pEsd->GetTrack(iTrk); Int_t ch=pTrk->GetHMPIDcluIdx(); //get track and chamber intersected by it | |
240 | if(ch<0) continue; //this track does not intersect any chamber | |
241 | Float_t thRa,phRa,xRa,yRa; pTrk->GetHMPIDtrk(xRa,yRa,thRa,phRa); //get info on current track | |
e4290ede | 242 | ch/=1000000; |
e9f11028 | 243 | Float_t xPc=0,yPc=0; AliHMPIDTracker::IntTrkCha(pTrk,xPc,yPc); //find again intersection of track with PC--> it is not stored in ESD! |
244 | fRenTxC[ch]->SetNextPoint(xPc,yPc); //add this intersection point | |
245 | Float_t ckov=pTrk->GetHMPIDsignal(); //get ckov angle stored for this track | |
45b05be5 | 246 | // Printf("RenderEsd: thetaC %f",ckov); |
e9f11028 | 247 | if(ckov>0){ |
248 | rec.SetTrack(xRa,yRa,thRa,phRa); | |
45b05be5 | 249 | for(Int_t j=0;j<100;j++){ |
e9f11028 | 250 | TVector2 pos; pos=rec.TracePhot(ckov,j*0.0628); |
ae5a42aa | 251 | if(!AliHMPIDParam::IsInDead(pos.X(),pos.Y())) fRenRin[ch]->SetNextPoint(pos.X(),pos.Y()); |
e9f11028 | 252 | } |
253 | }//if ckov is valid | |
254 | }//tracks loop | |
255 | }//RenEsd() | |
a1d55ff3 | 256 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
6e06db1d | 257 | void SimulateEsd(AliESDEvent *pEsd) |
a1d55ff3 | 258 | { |
259 | TParticle part; TLorentzVector mom; | |
50e69536 | 260 | for(Int_t iTrk=0;iTrk<100;iTrk++){//stack loop |
45b05be5 | 261 | Printf("SimulateEsd: kProton %i",kProton); |
a1d55ff3 | 262 | part.SetPdgCode(kProton); |
50e69536 | 263 | part.SetProductionVertex(0,0,0,0); |
264 | Double_t eta= -0.2+gRandom->Rndm()*0.4; //rapidity is random [-0.2,+0.2] | |
265 | Double_t phi= gRandom->Rndm()*60.*TMath::DegToRad(); //phi is random [ 0 , 60 ] degrees | |
266 | mom.SetPtEtaPhiM(5,eta,phi,part.GetMass()); | |
a1d55ff3 | 267 | part.SetMomentum(mom); |
268 | AliESDtrack trk(&part); | |
269 | pEsd->AddTrack(&trk); | |
270 | }//stack loop | |
271 | }//EsdFromStack() | |
272 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
6e06db1d | 273 | void SimulateHits(AliESDEvent *pEsd, TClonesArray *pHits) |
e9f11028 | 274 | {//used by SimulateEvent to simulate hits out from provided ESD |
50e69536 | 275 | const Int_t kCerenkov=50000050; |
276 | const Int_t kFeedback=50000051; | |
f3bae3e2 | 277 | |
a1d55ff3 | 278 | AliHMPIDRecon rec; |
e9f11028 | 279 | Float_t eMip=200e-9,ePho=7.5e-9; |
a591e55f | 280 | Int_t hc=0; |
a1d55ff3 | 281 | for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//tracks loop |
282 | AliESDtrack *pTrk=pEsd->GetTrack(iTrk); | |
e9f11028 | 283 | Float_t xRa,yRa; |
284 | Int_t ch=AliHMPIDTracker::IntTrkCha(pTrk,xRa,yRa); | |
a1d55ff3 | 285 | if(ch<0) continue; //this track does not hit HMPID |
50e69536 | 286 | Float_t beta = pTrk->GetP()/(TMath::Sqrt(pTrk->GetP()*pTrk->GetP()+0.938*0.938)); |
287 | Float_t ckov=TMath::ACos(1./(beta*1.292)); | |
a1d55ff3 | 288 | |
e9f11028 | 289 | Float_t theta,phi,xPc,yPc,; pTrk->GetHMPIDtrk(xPc,yPc,theta,phi); rec.SetTrack(xRa,yRa,theta,phi); |
a1d55ff3 | 290 | |
ae5a42aa | 291 | if(!AliHMPIDParam::IsInDead(xPc,yPc)) new((*pHits)[hc++]) AliHMPIDHit(ch,eMip,kProton ,iTrk,xPc,yPc); //mip hit |
50e69536 | 292 | Int_t nPhots = (Int_t)(20.*TMath::Power(TMath::Sin(ckov),2)/TMath::Power(TMath::Sin(TMath::ACos(1./1.292)),2)); |
293 | for(int i=0;i<nPhots;i++){ | |
a591e55f | 294 | TVector2 pos; |
295 | pos=rec.TracePhot(ckov,gRandom->Rndm()*TMath::TwoPi()); | |
ae5a42aa | 296 | if(!AliHMPIDParam::IsInDead(pos.X(),pos.Y())) new((*pHits)[hc++]) AliHMPIDHit(ch,ePho,kCerenkov,iTrk,pos.X(),pos.Y()); |
a1d55ff3 | 297 | } //photon hits |
e9f11028 | 298 | for(int i=0;i<3;i++){//feedback photons |
299 | Float_t x=gRandom->Rndm()*160; Float_t y=gRandom->Rndm()*150; | |
ae5a42aa | 300 | if(!AliHMPIDParam::IsInDead(x,y)) new((*pHits)[hc++]) AliHMPIDHit(ch,ePho,kFeedback,iTrk,x,y); //feedback hits |
e9f11028 | 301 | }//photon hits loop |
a1d55ff3 | 302 | }//tracks loop |
e9f11028 | 303 | }//SimulateHits() |
a1d55ff3 | 304 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
e9f11028 | 305 | void DoZoom(Int_t evt, Int_t px, Int_t py, TObject *) |
a1d55ff3 | 306 | { |
8282f9d0 | 307 | if(evt!=5 && evt!=6) return; //5- zoom in 6-zoom out |
308 | const Int_t minZoom=64; | |
309 | const Int_t maxZoom=2; | |
310 | static Int_t zoom=minZoom; //zoom level | |
311 | if(evt==5&&zoom==maxZoom) return; | |
312 | if(evt==6&&zoom==minZoom) return; | |
313 | ||
d73ece6a | 314 | // if(!obj->IsA()->InheritsFrom("TPad")) return; //current object is not pad |
50e69536 | 315 | TVirtualPad *pPad=gPad->GetSelectedPad(); |
8282f9d0 | 316 | if(pPad->GetNumber()==3 || pPad->GetNumber()==7) return; //current pad is wrong |
317 | ||
318 | // Printf("evt=%i (%i,%i) %s",evt,px,py,obj->GetName()); | |
319 | ||
320 | Float_t x=pPad->AbsPixeltoX(px); Float_t y=pPad->AbsPixeltoY(py); | |
6ed6a312 | 321 | |
322 | if(evt==5){ zoom=zoom/2; pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom in | |
323 | else { zoom=zoom*2; pPad->Range(x-zoom*2,y-zoom*2,x+zoom*2,y+zoom*2);} //zoom out | |
ae5a42aa | 324 | if(zoom==minZoom) pPad->Range(-10,-10,AliHMPIDParam::SizeAllX()+5,AliHMPIDParam::SizeAllY()+5); |
8282f9d0 | 325 | ((TCanvas *)gTQSender)->SetTitle(Form("zoom x%i",minZoom/zoom)); |
326 | pPad->Modified(); | |
327 | pPad->Update(); | |
a1d55ff3 | 328 | } |
8282f9d0 | 329 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
e9f11028 | 330 | void ReadEvent() |
331 | {//used by NextEvent() to read curent event and construct all render elements | |
332 | if(fNevt && fEvt>=fNevt) fEvt=0; //loop over max event | |
45b05be5 | 333 | if(fNevt && fEvt<0) fEvt=fNevt-1; //loop over max event |
d73ece6a | 334 | |
e9f11028 | 335 | if(fHitFile){ |
336 | if(fHitTree) delete fHitTree; fHitTree=(TTree*)fHitFile->Get(Form("Event%i/TreeH",fEvt)); | |
337 | ||
338 | fHitTree->SetBranchAddress("HMPID",&fHitLst); | |
339 | for(Int_t iEnt=0;iEnt<fHitTree->GetEntries();iEnt++){ | |
340 | fHitTree->GetEntry(iEnt); | |
341 | RenderHit(fHitLst); | |
342 | }//prim loop | |
343 | }//if hits file | |
e4290ede | 344 | |
345 | if(fDigFile){ | |
346 | if(fDigTree) delete fDigTree; fDigTree=(TTree*)fDigFile->Get(Form("Event%i/TreeD",fEvt)); | |
347 | for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) fDigTree->SetBranchAddress(Form("HMPID%i",iCh),&(*fDigLst)[iCh]); | |
348 | fDigTree->GetEntry(0); | |
349 | RenderDig(fDigLst); | |
350 | }//if digs file | |
351 | ||
e9f11028 | 352 | if(fCluFile){ |
353 | if(fCluTree) delete fCluTree; fCluTree=(TTree*)fCluFile->Get(Form("Event%i/TreeR",fEvt)); | |
ae5a42aa | 354 | for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) fCluTree->SetBranchAddress(Form("HMPID%i",iCh),&(*fCluLst)[iCh]); |
e9f11028 | 355 | fCluTree->GetEntry(0); |
e4290ede | 356 | RenderClu(fCluLst); |
e9f11028 | 357 | }//if clus file |
358 | ||
359 | if(fEsdFile){//if ESD file open | |
360 | fEsdTree->GetEntry(fEvt); | |
361 | RenderEsd(fEsd); | |
362 | }//if ESD file | |
363 | }//ReadEvent() | |
d73ece6a | 364 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
e9f11028 | 365 | void SimulateEvent() |
366 | {// used by NextEvent() to simulate all info | |
367 | AliCDBManager* pCDB = AliCDBManager::Instance(); pCDB->SetDefaultStorage("local://$HOME"); pCDB->SetRun(0); | |
368 | AliCDBEntry *pNmeanEnt=pCDB->Get("HMPID/Calib/Nmean"); | |
369 | ||
370 | SimulateEsd(fEsd); | |
371 | SimulateHits(fEsd,fHitLst); | |
45b05be5 | 372 | AliHMPIDv2::Hit2Sdi(fHitLst,fSdiLst); |
e9f11028 | 373 | AliHMPIDDigitizer::Sdi2Dig(fSdiLst,fDigLst); |
374 | AliHMPIDReconstructor::Dig2Clu(fDigLst,fCluLst); | |
375 | AliHMPIDTracker::Recon(fEsd,fCluLst,(TObjArray*)pNmeanEnt->GetObject()); | |
376 | ||
377 | RenderHit(fHitLst); | |
e4290ede | 378 | RenderClu(fCluLst); |
e9f11028 | 379 | RenderEsd(fEsd); |
380 | }//SimulateEvent() | |
d73ece6a | 381 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
45b05be5 | 382 | void GetEvent() |
d73ece6a | 383 | { |
e9f11028 | 384 | ClearRenders(); |
385 | switch(fType){ | |
386 | case 1: ReadEvent();break; | |
387 | // case 2: ReadCosmic(); break; | |
388 | case 3: SimulateEvent(); break; | |
389 | default: return; | |
390 | } | |
391 | Draw(); | |
d73ece6a | 392 | } |
393 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
45b05be5 | 394 | void NextEvent() |
395 | { | |
396 | fTimeArrow = 1; | |
397 | fEvt+=fTimeArrow; | |
398 | GetEvent(); | |
399 | } | |
400 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
401 | void PrevEvent() | |
402 | { | |
403 | fTimeArrow = -1; | |
404 | fEvt+=fTimeArrow; | |
405 | GetEvent(); | |
406 | } | |
407 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
408 | void End() | |
409 | { | |
410 | gSystem->Exit(0); | |
411 | } | |
412 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
e9f11028 | 413 | void Hdisp() |
414 | {//display events from files if any in current directory or simulated events | |
b6eacec7 | 415 | AliHMPIDParam::Instance(); // first invocation of AliHMPIDParam to initialize geometry... |
e9f11028 | 416 | CreateContainers(); |
417 | CreateRenders(); | |
d73ece6a | 418 | |
e9f11028 | 419 | TString title="Session with"; |
420 | if(gSystem->IsFileInIncludePath("HMPID.Hits.root")){// tries to open hits | |
421 | fHitFile=TFile::Open("HMPID.Hits.root"); fNevt=fHitFile->GetNkeys(); fType=1; title+=Form(" HITS-%i ",fNevt); | |
422 | } | |
6b155609 | 423 | |
e9f11028 | 424 | if(gSystem->IsFileInIncludePath("HMPID.Digits.root")){// tries to open clusters |
e4290ede | 425 | fDigFile=TFile::Open("HMPID.Digits.root"); fNevt=fDigFile->GetNkeys(); fType=1; title+=Form(" DIGITS-%i ",fNevt); |
426 | } | |
6b155609 | 427 | |
e9f11028 | 428 | if(gSystem->IsFileInIncludePath("HMPID.RecPoints.root")){// tries to open clusters |
429 | fCluFile=TFile::Open("HMPID.RecPoints.root"); fNevt=fCluFile->GetNkeys(); fType=1; title+=Form(" CLUSTERS-%i ",fNevt); | |
e4290ede | 430 | } |
50e69536 | 431 | |
e9f11028 | 432 | if(gSystem->IsFileInIncludePath("AliESDs.root")){ |
e4290ede | 433 | fEsdFile=TFile::Open("AliESDs.root"); |
434 | fEsdTree=(TTree*)fEsdFile->Get("esdTree"); | |
435 | fEsd->ReadFromTree(fEsdTree); fEsd->GetStdContent(); //clm: new ESD schema: see Task Force meeting 20th June, 2007 | |
436 | fNevt=fEsdTree->GetEntries(); fType=1; title+=Form(" ESD-%i ",fNevt); | |
437 | ||
438 | //clm: we need to set the magnetic field | |
439 | if(gSystem->IsFileInIncludePath("galice.root")){ | |
440 | if(gAlice) delete gAlice; | |
441 | gAL=AliRunLoader::Open(); | |
442 | gAL->LoadgAlice(); | |
443 | AliHMPIDTracker::SetFieldMap(gAL->GetAliRun()->Field(),kTRUE); | |
444 | }else{ | |
445 | Printf("=============== NO galice file! Magnetic field for ESD tracking is: %f ===============",AliTracker::GetBz()); | |
446 | } | |
447 | ||
e9f11028 | 448 | } |
449 | ||
450 | if(gSystem->IsFileInIncludePath("cosmic.root")){ //clm: Check if cosmic file is in the folder | |
451 | fCosFile=TFile::Open("cosmic.root"); fCosTree=(TTree*)fCosFile->Get("cosmic"); fNevt=fCosTree->GetEntries(); fType=2; | |
452 | fCosTree->SetBranchAddress("Digs",&fDigLst); fCosTree->SetBranchAddress("Clus",&fCluLst); | |
453 | } | |
454 | ||
45b05be5 | 455 | fCanvas=new TCanvas("all","",-1024,768);fCanvas->SetWindowSize(1024,768); |
456 | fCanvas->Divide(3,3,0,0);// pAll->ToggleEditor(); | |
e9f11028 | 457 | fCanvas->Connect("ProcessedEvent(Int_t,Int_t,Int_t,TObject*)",0,0,"DoZoom(Int_t,Int_t,Int_t,TObject*)"); |
458 | fCanvas->cd(7); | |
459 | switch(fType){ | |
460 | case 1: fCanvas->SetTitle(title.Data()); | |
45b05be5 | 461 | TButton *pNxtBtn = new TButton("Next" ,"NextEvent()" ,0.0,0.0,0.3,0.1); pNxtBtn->Draw(); |
462 | TButton *pPrvBtn = new TButton("Previous" ,"PrevEvent()" ,0.3,0.0,0.6,0.1); pPrvBtn->Draw(); | |
463 | if(fHitFile){TButton *pHitBtn = new TButton("Print hits","PrintHits()" ,0.0,0.2,0.3,0.3); pHitBtn->Draw();} | |
464 | if(fDigFile){TButton *pDigBtn = new TButton("Print digs","PrintDigs()" ,0.0,0.4,0.3,0.5); pDigBtn->Draw();} | |
465 | if(fCluFile){TButton *pCluBtn = new TButton("Print clus","PrintClus()" ,0.0,0.6,0.3,0.7); pCluBtn->Draw();} | |
466 | if(fEsdFile){TButton *pEsdBtn = new TButton("Print ESD ","PrintEsd()" ,0.0,0.8,0.3,0.9); pEsdBtn->Draw();} | |
467 | if(fHitFile){ fHitMipBok = new TButton(Form("Mip %s",fStHitMip.Data()),"SwitchHitMip()",0.3,0.16,0.6,0.22); fHitMipBok->Draw();} | |
468 | if(fHitFile){ fHitCkoBok = new TButton(Form("Ckov %s",fStHitCko.Data()),"SwitchHitCko()",0.3,0.22,0.6,0.28); fHitCkoBok->Draw();} | |
469 | if(fHitFile){ fHitFeeBok = new TButton(Form("Fdbk %s",fStHitFee.Data()),"SwitchHitFee()",0.3,0.28,0.6,0.34); fHitFeeBok->Draw();} | |
470 | if(fDigFile){ fDigBok = new TButton(fStDig ,"SwitchDigs()" ,0.3,0.4,0.6,0.5); fDigBok->Draw();} | |
471 | if(fCluFile){ fCluBok = new TButton(fStClu ,"SwitchClus()" ,0.3,0.6,0.6,0.7); fCluBok->Draw();} | |
472 | if(fEsdFile){ fEsdBok = new TButton(fStEsd ,"SwitchEsd()" ,0.3,0.8,0.6,0.9); fEsdBok->Draw();} | |
473 | TButton *pEnd = new TButton("Quit" ,"End()" ,0.6,0.0,0.9,0.1); pEnd->Draw(); | |
e9f11028 | 474 | break; |
475 | case 2: fCanvas->SetTitle("COSMIC"); | |
45b05be5 | 476 | TButton *pCosBtn = new TButton("Next" ,"NextEvent()",0.0,0.0,0.3,0.1); pCosBtn->Draw(); |
e9f11028 | 477 | |
478 | break; | |
479 | case 3: fCanvas->SetTitle("SIMULATION"); | |
45b05be5 | 480 | TButton *pSimBtn = new TButton("Simulate" ,"NextEvent()",0.0,0.0,0.3,0.1); pSimBtn->Draw(); |
481 | TButton *pHitBtn = new TButton("Print hits","PrintHits()",0.0,0.2,0.3,0.3); pHitBtn->Draw(); | |
482 | TButton *pDigBtn = new TButton("Print digs","PrintDigs()",0.0,0.4,0.3,0.5); pDigBtn->Draw(); | |
483 | TButton *pCluBtn = new TButton("Print clus","PrintClus()",0.0,0.6,0.3,0.7); pCluBtn->Draw(); | |
484 | TButton *pEsdBtn = new TButton("Print ESD" ,"PrintEsd()" ,0.0,0.8,0.3,0.9); pEsdBtn->Draw(); | |
485 | if(fHitFile){ fHitMipBok = new TButton(Form("Mips %s",fStHitMip.Data()),"SwitchHitMip()",0.3,0.16,0.6,0.22); fHitMipBok->Draw();} | |
486 | if(fHitFile){ fHitCkoBok = new TButton(Form("Ckov %s",fStHitCko.Data()),"SwitchHitCko()",0.3,0.22,0.6,0.28); fHitCkoBok->Draw();} | |
487 | if(fHitFile){ fHitFeeBok = new TButton(Form("Fdbk %s",fStHitFee.Data()),"SwitchHitFee()",0.3,0.28,0.6,0.34); fHitFeeBok->Draw();} | |
488 | if(fDigFile){ fDigBok = new TButton(fStDig ,"SwitchDigs()" ,0.3,0.4,0.6,0.5); fDigBok->Draw();} | |
489 | if(fCluFile){ fCluBok = new TButton(fStClu ,"SwitchClus()" ,0.3,0.6,0.6,0.7); fCluBok->Draw();} | |
490 | if(fEsdFile){ fEsdBok = new TButton(fStEsd ,"SwitchEsd()" ,0.3,0.8,0.6,0.9); fEsdBok->Draw();} | |
e9f11028 | 491 | break; |
492 | } | |
493 | ||
0112d781 | 494 | NextEvent(); |
e9f11028 | 495 | |
496 | }//Hdisp() | |
497 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
45b05be5 | 498 | void SwitchHitMip() |
499 | { | |
500 | fStHitMip=(fStHitMip=="OFF")? "ON":"OFF"; | |
501 | GetEvent(); | |
502 | } | |
503 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
504 | void SwitchHitCko() | |
505 | { | |
506 | fStHitCko=(fStHitCko=="OFF")? "ON":"OFF"; | |
507 | GetEvent(); | |
508 | } | |
509 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
510 | void SwitchHitFee() | |
511 | { | |
512 | fStHitFee=(fStHitFee=="OFF")? "ON":"OFF"; | |
513 | GetEvent(); | |
514 | } | |
515 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
516 | void SwitchDigs() | |
517 | { | |
518 | fStDig=(fStDig=="OFF")? "ON":"OFF"; | |
519 | GetEvent(); | |
520 | } | |
521 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
522 | void SwitchClus() | |
523 | { | |
524 | fStClu=(fStClu=="OFF")? "ON":"OFF"; | |
525 | GetEvent(); | |
526 | } | |
527 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ | |
528 | void SwitchEsd() | |
529 | { | |
530 | fStEsd=(fStEsd=="OFF")? "ON":"OFF"; | |
531 | GetEvent(); | |
532 | } | |
533 | //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |