1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 Revision 1.5 2001/01/26 19:58:48 hristov
19 Major upgrade of AliRoot code
21 Revision 1.4 2001/01/17 10:50:50 hristov
22 Corrections to destructors
24 Revision 1.3 2000/12/18 14:16:31 alibrary
27 Revision 1.2 2000/12/18 11:33:48 alibrary
28 New call frequence histograms per module and volume
30 Revision 1.1 2000/11/30 07:12:48 alibrary
31 Introducing new Rndm and QA classes
35 ///////////////////////////////////////////////////////////////////////////////
38 ///////////////////////////////////////////////////////////////////////////////
43 #include "TObjArray.h"
49 #include "TLorentzVector.h"
50 #include "TDatabasePDG.h"
54 #include "TPaveLabel.h"
59 #include "AliModule.h"
65 //_____________________________________________________________________________
66 AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
67 fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
70 // Default constructor
74 //_____________________________________________________________________________
75 AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
78 // Constructor, creates the list of lists of histograms
86 fQAList = new TObjArray(ndets);
87 TObjArray &hist = *fQAList;
91 TObjArray &mods = *(gAlice->Modules());
93 TList *dir = gDirectory->GetList();
94 for (i=0; i<ndets; i++) {
95 hist[i] = list = new TList();
96 mod = (AliModule *) mods[i];
99 sprintf(title,"Spectrum entering: %s ",mod->GetName());
100 list->Add(new TH1F("hEnIn",strcat(title,mod->GetTitle()),100,-4,2));
101 dir->Remove(dir->FindObject("hEnIn"));
103 sprintf(title,"Spectrum exiting: %s ",mod->GetName());
104 list->Add(new TH1F("hEnOut",strcat(title,mod->GetTitle()),100,-4,2));
105 dir->Remove(dir->FindObject("hEnOut"));
108 sprintf(title,"Z coordinate entering: %s ",mod->GetName());
109 list->Add(new TH1F("hZIn",strcat(title,mod->GetTitle()),100,
110 mod->ZMin(),mod->ZMax()));
111 dir->Remove(dir->FindObject("hZIn"));
113 sprintf(title,"Z coordinate exiting: %s ",mod->GetName());
114 list->Add(new TH1F("hZOut",strcat(title,mod->GetTitle()),100,
115 mod->ZMin(),mod->ZMax()));
116 dir->Remove(dir->FindObject("hZOut"));
119 gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
121 fDetDone = new Int_t[fNdets];
124 // Global QA histograms
126 fQAHist = new TObjArray(2);
127 fNvolumes=gMC->NofVolumes();
129 fQAHist->Add(new TH1F("hMCVcalls","Monte Carlo calls per volume",
130 fNvolumes, 0.5, fNvolumes+0.5));
131 h = (TH1F*) dir->FindObject("hMCVcalls");
132 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddVolumeName()"));
133 dir->Remove(dir->FindObject("hMCVcalls"));
135 // Build list of volume names
137 fVolNames=new TObjArray(fNvolumes);
138 for(i=0;i<fNvolumes;++i) {
139 AliModule *mod = (AliModule*)
140 (*gAlice->Modules())[gAlice->DetFromMate(gMC->VolId2Mate(i+1))];
141 (*fVolNames)[i]=new TNamed(gMC->VolName(i+1),mod->GetName());
144 fQAHist->Add(new TH1F("hMCMcalls","Monte Carlo calls per module",
145 fNdets, -0.5, fNdets-0.5));
146 h = (TH1F*) dir->FindObject("hMCMcalls");
147 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddModuleName()"));
149 dir->Remove(dir->FindObject("hMCMcalls"));
151 // Build list of module names
153 fModNames=new TObjArray(fNdets);
154 for(i=0;i<fNdets;++i)
156 new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
159 //_____________________________________________________________________________
161 AliMCQA::~AliMCQA() {
187 //_____________________________________________________________________________
188 void AliMCQA::Browse(TBrowser *b)
191 // Called when the item "Run" is clicked on the left pane
192 // of the Root browser.
193 // It displays the Root Trees and all detectors.
197 // Global histos first
199 TIter global(fQAHist);
200 while((hist = (TH1*)global()))
201 b->Add(hist,hist->GetTitle());
203 // Module histograms now
207 while((histos = (TList*)next())) {
209 while((hist = (TH1*)next1()))
210 b->Add(hist,hist->GetTitle());
214 //_____________________________________________________________________________
215 void AliMCQA::PreTrack()
218 // Called before each track
221 for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
224 //_____________________________________________________________________________
225 void AliMCQA::StepManager(Int_t id)
228 // Called at each step
233 // Fill Global histograms first
237 static TH1F* mcvcalls = (TH1F*) fQAHist->FindObject("hMCVcalls");
238 mcvcalls->Fill(gMC->CurrentVolID(copy));
239 static TH1F* mcmcalls = (TH1F*) fQAHist->FindObject("hMCMcalls");
243 // Now the step manager histograms
246 static Double_t mpi0=0;
247 static Double_t mpip=0;
248 static Double_t mpim=0;
249 static Double_t mep=0;
250 static Double_t mem=0;
252 Int_t num = gMC->TrackPid();
256 if (mpi0==0) mpi0=gAlice->PDGDB()->GetParticle(num)->Mass();
260 if (mpip==0) mpip=gAlice->PDGDB()->GetParticle(num)->Mass();
264 if (mpim==0) mpim=gAlice->PDGDB()->GetParticle(num)->Mass();
268 if (mep==0) mep=gAlice->PDGDB()->GetParticle(num)->Mass();
272 if (mem==0) mem=gAlice->PDGDB()->GetParticle(num)->Mass();
276 mass =gAlice->PDGDB()->GetParticle(num)->Mass();
280 static TLorentzVector p, x;
281 gMC->TrackMomentum(p);
282 gMC->TrackPosition(x);
283 Double_t energy = TMath::Max(p[3]-mass,1.e-12);
285 if(!fDetDone[fOldId] && !gMC->IsNewTrack()) {
286 TList *histold = (TList*) (*fQAList)[fOldId];
287 hist = (TH1F*) histold->FindObject("hEnOut");
288 hist->Fill(TMath::Log10(energy));
289 hist = (TH1F*) histold->FindObject("hZOut");
294 if(!fDetDone[id] && !gMC->IsNewTrack()) {
295 TList *histnew = (TList*) (*fQAList)[id];
296 hist = (TH1F*) histnew->FindObject("hEnIn");
297 hist->Fill(TMath::Log10(energy));
298 hist = (TH1F*) histnew->FindObject("hZIn");
305 //_____________________________________________________________________________
306 void AliMCQA::AddModuleName()
309 // Add function DrawModuleName to the module frequency histogram pad
311 static TVirtualPad* oldpad=0;
313 gPad->GetCanvas()->FeedbackMode(kTRUE);
314 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
315 DrawPaveLabel(fMPaveLabel);
320 //_____________________________________________________________________________
321 void AliMCQA::AddVolumeName()
324 // Add function DrawVolumeName to the volume frequency histogram pad
326 static TVirtualPad* oldpad=0;
328 gPad->GetCanvas()->FeedbackMode(kTRUE);
329 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
330 DrawPaveLabel(fVPaveLabel);
335 //_____________________________________________________________________________
336 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
339 // Draws the PaveLabel with the meaning of the bin
341 float uxmin = gPad->GetUxmin();
342 float uxmax = gPad->GetUxmax();
343 float uymin = gPad->GetUymin();
344 float uymax = gPad->GetUymax();
345 float lx = uxmax-uxmin;
346 float ly = uymax-uymin;
350 = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
351 uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
356 //_____________________________________________________________________________
357 Int_t AliMCQA::GetHBin(const char* hname)
360 // Get the bin where the cursor is
362 TList *dir = gDirectory->GetList();
363 TH1 *h=(TH1*)dir->FindObject(hname);
366 int px = gPad->GetEventX();
367 Float_t upx = gPad->AbsPixeltoX(px);
368 Float_t x = gPad->PadtoX(upx);
370 return h->GetXaxis()->FindBin(x);
373 //_____________________________________________________________________________
374 void AliMCQA::DrawModuleName()
377 // Writes the name of the module of the bin where the cursor is
379 TObject *select = gPad->GetSelected();
382 Int_t binx = GetHBin("hMCMcalls");
383 if(0<binx && binx<=fNdets) {
385 strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
386 fMPaveLabel->SetLabel(lab);
393 //_____________________________________________________________________________
394 void AliMCQA::DrawVolumeName()
397 // Writes the name of the volume:module of the bin where the cursor is
399 TObject *select = gPad->GetSelected();
402 Int_t binx = GetHBin("hMCVcalls");
403 if(0<binx && binx<=fNvolumes) {
405 sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
406 ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
407 fVPaveLabel->SetLabel(lab);