]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliMCQA.cxx
Fix for multiple events per file: inhibit decrease of size of fParticleFileMap.
[u/mrichter/AliRoot.git] / STEER / AliMCQA.cxx
CommitLineData
65fb704d 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
16/*
17$Log$
2ab0c725 18Revision 1.4 2001/01/17 10:50:50 hristov
19Corrections to destructors
20
e460afec 21Revision 1.3 2000/12/18 14:16:31 alibrary
22HP compatibility fix
23
af2e514f 24Revision 1.2 2000/12/18 11:33:48 alibrary
25New call frequence histograms per module and volume
26
8409d8c9 27Revision 1.1 2000/11/30 07:12:48 alibrary
28Introducing new Rndm and QA classes
29
65fb704d 30*/
31
32///////////////////////////////////////////////////////////////////////////////
33// //
34// //
35///////////////////////////////////////////////////////////////////////////////
36
8409d8c9 37
65fb704d 38#include <strings.h>
39
40#include "TObjArray.h"
41#include "TH1.h"
42#include "TList.h"
43#include "TROOT.h"
44#include "TBrowser.h"
45#include "TMath.h"
46#include "TLorentzVector.h"
47#include "TDatabasePDG.h"
48#include "TMath.h"
8409d8c9 49#include "TPad.h"
50#include "TExec.h"
51#include "TPaveLabel.h"
52#include "TCanvas.h"
65fb704d 53
54#include "AliMCQA.h"
55#include "AliRun.h"
56#include "AliModule.h"
57#include "AliMC.h"
58
59ClassImp(AliMCQA)
60
61
62//_____________________________________________________________________________
8409d8c9 63AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
64 fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
65{
66 //
67 // Default constructor
68 //
69}
70
71//_____________________________________________________________________________
72AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
65fb704d 73{
74 //
75 // Constructor, creates the list of lists of histograms
76 //
77 TList *list;
8409d8c9 78 TH1F* h;
af2e514f 79 Int_t i;
8409d8c9 80
81 fNdets=ndets;
65fb704d 82
83 fQAList = new TObjArray(ndets);
84 TObjArray &hist = *fQAList;
85
86 char title[100];
87 //
88 TObjArray &mods = *(gAlice->Modules());
89 AliModule *mod;
90 TList *dir = gDirectory->GetList();
af2e514f 91 for (i=0; i<ndets; i++) {
65fb704d 92 hist[i] = list = new TList();
93 mod = (AliModule *) mods[i];
94
95 // Energy Spectrum
96 sprintf(title,"Spectrum entering: %s ",mod->GetName());
97 list->Add(new TH1F("hEnIn",strcat(title,mod->GetTitle()),100,-4,2));
98 dir->Remove(dir->FindObject("hEnIn"));
99
100 sprintf(title,"Spectrum exiting: %s ",mod->GetName());
101 list->Add(new TH1F("hEnOut",strcat(title,mod->GetTitle()),100,-4,2));
102 dir->Remove(dir->FindObject("hEnOut"));
103
104 // Z position
105 sprintf(title,"Z coordinate entering: %s ",mod->GetName());
106 list->Add(new TH1F("hZIn",strcat(title,mod->GetTitle()),100,
107 mod->ZMin(),mod->ZMax()));
108 dir->Remove(dir->FindObject("hZIn"));
109
110 sprintf(title,"Z coordinate exiting: %s ",mod->GetName());
111 list->Add(new TH1F("hZOut",strcat(title,mod->GetTitle()),100,
112 mod->ZMin(),mod->ZMax()));
113 dir->Remove(dir->FindObject("hZOut"));
114 }
115 //
116 gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
117
118 fDetDone = new Int_t[fNdets];
8409d8c9 119
120 //
121 // Global QA histograms
122 //
123 fQAHist = new TObjArray(2);
124 fNvolumes=gMC->NofVolumes();
125
126 fQAHist->Add(new TH1F("hMCVcalls","Monte Carlo calls per volume",
127 fNvolumes, 0.5, fNvolumes+0.5));
128 h = (TH1F*) dir->FindObject("hMCVcalls");
129 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddVolumeName()"));
130 dir->Remove(dir->FindObject("hMCVcalls"));
131 //
132 // Build list of volume names
133 //
8409d8c9 134 fVolNames=new TObjArray(fNvolumes);
135 for(i=0;i<fNvolumes;++i) {
136 AliModule *mod = (AliModule*)
137 (*gAlice->Modules())[gAlice->DetFromMate(gMC->VolId2Mate(i+1))];
138 (*fVolNames)[i]=new TNamed(gMC->VolName(i+1),mod->GetName());
139 }
140
141 fQAHist->Add(new TH1F("hMCMcalls","Monte Carlo calls per module",
142 fNdets, -0.5, fNdets-0.5));
143 h = (TH1F*) dir->FindObject("hMCMcalls");
144 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddModuleName()"));
145
146 dir->Remove(dir->FindObject("hMCMcalls"));
147 //
148 // Build list of module names
149 //
150 fModNames=new TObjArray(fNdets);
151 for(i=0;i<fNdets;++i)
152 (*fModNames)[i]=
153 new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
65fb704d 154}
155
e460afec 156//_____________________________________________________________________________
157
158AliMCQA::~AliMCQA() {
159 if (fQAList) {
160 fQAList->Delete();
161 delete fQAList;
162 fQAList=0;
163 }
164 if (fQAHist) {
165 fQAHist->Delete();
166 delete fQAHist;
167 fQAHist=0;
168 }
169 if (fVolNames) {
170 fVolNames->Delete();
171 delete fVolNames;
172 fVolNames=0;
173 }
174 if (fModNames) {
175 fModNames->Delete();
176 delete fModNames;
177 fModNames=0;
178 }
179}
180
65fb704d 181//_____________________________________________________________________________
182void AliMCQA::Browse(TBrowser *b)
183{
184 //
185 // Called when the item "Run" is clicked on the left pane
186 // of the Root browser.
187 // It displays the Root Trees and all detectors.
188 //
8409d8c9 189 TH1 *hist;
190 //
191 // Global histos first
192 //
193 TIter global(fQAHist);
194 while((hist = (TH1*)global()))
195 b->Add(hist,hist->GetTitle());
196 //
197 // Module histograms now
198 //
65fb704d 199 TIter next(fQAList);
200 TList *histos;
65fb704d 201 while((histos = (TList*)next())) {
202 TIter next1(histos);
8409d8c9 203 while((hist = (TH1*)next1()))
65fb704d 204 b->Add(hist,hist->GetTitle());
65fb704d 205 }
206}
207
208//_____________________________________________________________________________
209void AliMCQA::PreTrack()
210{
8409d8c9 211 //
212 // Called before each track
213 //
65fb704d 214 fOldId=-1;
215 for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
216}
217
218//_____________________________________________________________________________
219void AliMCQA::StepManager(Int_t id)
220{
8409d8c9 221 //
222 // Called at each step
223 //
224 TH1F *hist;
225 Int_t copy;
226 //
227 // Fill Global histograms first
228 //
2ab0c725 229
230
231 static TH1F* mcvcalls = (TH1F*) fQAHist->FindObject("hMCVcalls");
232 mcvcalls->Fill(gMC->CurrentVolID(copy));
233 static TH1F* mcmcalls = (TH1F*) fQAHist->FindObject("hMCMcalls");
234 mcmcalls->Fill(id);
235
8409d8c9 236 //
237 // Now the step manager histograms
238 //
65fb704d 239 if(fOldId != id) {
2ab0c725 240 static Double_t mpi0=0;
241 static Double_t mpip=0;
242 static Double_t mpim=0;
243 static Double_t mep=0;
244 static Double_t mem=0;
245 Double_t mass=0;
246 Int_t num = gMC->TrackPid();
247
248 switch (num) {
249 case 111:
250 if (mpi0==0) mpi0=gAlice->PDGDB()->GetParticle(num)->Mass();
251 mass=mpi0;
252 break;
253 case 211:
254 if (mpip==0) mpip=gAlice->PDGDB()->GetParticle(num)->Mass();
255 mass=mpip;
256 break;
257 case -211:
258 if (mpim==0) mpim=gAlice->PDGDB()->GetParticle(num)->Mass();
259 mass=mpim;
260 break;
261 case 11:
262 if (mep==0) mep=gAlice->PDGDB()->GetParticle(num)->Mass();
263 mass=mep;
264 break;
265 case -11:
266 if (mem==0) mem=gAlice->PDGDB()->GetParticle(num)->Mass();
267 mass=mem;
268 break;
269 default:
270 mass =gAlice->PDGDB()->GetParticle(num)->Mass();
271 break;
272 }
273
274 static TLorentzVector p, x;
65fb704d 275 gMC->TrackMomentum(p);
276 gMC->TrackPosition(x);
2ab0c725 277 Double_t energy = TMath::Max(p[3]-mass,1.e-12);
65fb704d 278 if(fOldId > -1) {
279 if(!fDetDone[fOldId] && !gMC->IsNewTrack()) {
280 TList *histold = (TList*) (*fQAList)[fOldId];
281 hist = (TH1F*) histold->FindObject("hEnOut");
282 hist->Fill(TMath::Log10(energy));
283 hist = (TH1F*) histold->FindObject("hZOut");
284 hist->Fill(x[2]);
285 fDetDone[fOldId]=1;
286 }
287 }
288 if(!fDetDone[id] && !gMC->IsNewTrack()) {
289 TList *histnew = (TList*) (*fQAList)[id];
290 hist = (TH1F*) histnew->FindObject("hEnIn");
291 hist->Fill(TMath::Log10(energy));
292 hist = (TH1F*) histnew->FindObject("hZIn");
293 hist->Fill(x[2]);
294 }
295 fOldId=id;
296 }
297}
8409d8c9 298
299//_____________________________________________________________________________
300void AliMCQA::AddModuleName()
301{
302 //
303 // Add function DrawModuleName to the module frequency histogram pad
304 //
305 static TVirtualPad* oldpad=0;
306 if(oldpad!=gPad) {
307 gPad->GetCanvas()->FeedbackMode(kTRUE);
308 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
309 DrawPaveLabel(fMPaveLabel);
310 oldpad=gPad;
311 }
312}
313
314//_____________________________________________________________________________
315void AliMCQA::AddVolumeName()
316{
317 //
318 // Add function DrawVolumeName to the volume frequency histogram pad
319 //
320 static TVirtualPad* oldpad=0;
321 if(oldpad!=gPad) {
322 gPad->GetCanvas()->FeedbackMode(kTRUE);
323 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
324 DrawPaveLabel(fVPaveLabel);
325 oldpad=gPad;
326 }
327}
328
329//_____________________________________________________________________________
330void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
331{
332 //
333 // Draws the PaveLabel with the meaning of the bin
334 //
335 float uxmin = gPad->GetUxmin();
336 float uxmax = gPad->GetUxmax();
337 float uymin = gPad->GetUymin();
338 float uymax = gPad->GetUymax();
339 float lx = uxmax-uxmin;
340 float ly = uymax-uymin;
341
342 if(pv) delete pv;
343 pv
344 = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
345 uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
346 "");
347 pv->Draw();
348}
349
350//_____________________________________________________________________________
351Int_t AliMCQA::GetHBin(const char* hname)
352{
353 //
354 // Get the bin where the cursor is
355 //
356 TList *dir = gDirectory->GetList();
357 TH1 *h=(TH1*)dir->FindObject(hname);
358
359
360 int px = gPad->GetEventX();
361 Float_t upx = gPad->AbsPixeltoX(px);
362 Float_t x = gPad->PadtoX(upx);
363
364 return h->GetXaxis()->FindBin(x);
365}
366
367//_____________________________________________________________________________
368void AliMCQA::DrawModuleName()
369{
370 //
371 // Writes the name of the module of the bin where the cursor is
372 //
373 TObject *select = gPad->GetSelected();
374 if(!select) return;
375
376 Int_t binx = GetHBin("hMCMcalls");
377 if(0<binx && binx<=fNdets) {
378 char lab[15];
379 strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
380 fMPaveLabel->SetLabel(lab);
381
382 gPad->Modified();
383 gPad->Update();
384 }
385}
386
387//_____________________________________________________________________________
388void AliMCQA::DrawVolumeName()
389{
390 //
391 // Writes the name of the volume:module of the bin where the cursor is
392 //
393 TObject *select = gPad->GetSelected();
394 if(!select) return;
395
396 Int_t binx = GetHBin("hMCVcalls");
397 if(0<binx && binx<=fNvolumes) {
398 char lab[20];
399 sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
400 ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
401 fVPaveLabel->SetLabel(lab);
402
403 gPad->Modified();
404 gPad->Update();
405 }
406}
407
408