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.6.6.1 2002/10/12 21:41:00 hristov
19 Remove the object from the list of browsables
21 Revision 1.6 2001/12/05 08:31:25 hristov
22 Destructor corrected, thanks to R.Brun
24 Revision 1.5 2001/01/26 19:58:48 hristov
25 Major upgrade of AliRoot code
27 Revision 1.4 2001/01/17 10:50:50 hristov
28 Corrections to destructors
30 Revision 1.3 2000/12/18 14:16:31 alibrary
33 Revision 1.2 2000/12/18 11:33:48 alibrary
34 New call frequence histograms per module and volume
36 Revision 1.1 2000/11/30 07:12:48 alibrary
37 Introducing new Rndm and QA classes
41 ///////////////////////////////////////////////////////////////////////////////
44 ///////////////////////////////////////////////////////////////////////////////
49 #include "TObjArray.h"
55 #include "TLorentzVector.h"
56 #include "TDatabasePDG.h"
60 #include "TPaveLabel.h"
65 #include "AliModule.h"
71 //_____________________________________________________________________________
72 AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
73 fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
76 // Default constructor
80 //_____________________________________________________________________________
81 AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
84 // Constructor, creates the list of lists of histograms
92 fQAList = new TObjArray(ndets);
93 TObjArray &hist = *fQAList;
97 TObjArray &mods = *(gAlice->Modules());
99 TList *dir = gDirectory->GetList();
100 for (i=0; i<ndets; i++) {
101 hist[i] = list = new TList();
102 mod = (AliModule *) mods[i];
105 sprintf(title,"Spectrum entering: %s ",mod->GetName());
106 list->Add(new TH1F("hEnIn",strcat(title,mod->GetTitle()),100,-4,2));
107 dir->Remove(dir->FindObject("hEnIn"));
109 sprintf(title,"Spectrum exiting: %s ",mod->GetName());
110 list->Add(new TH1F("hEnOut",strcat(title,mod->GetTitle()),100,-4,2));
111 dir->Remove(dir->FindObject("hEnOut"));
114 sprintf(title,"Z coordinate entering: %s ",mod->GetName());
115 list->Add(new TH1F("hZIn",strcat(title,mod->GetTitle()),100,
116 mod->ZMin(),mod->ZMax()));
117 dir->Remove(dir->FindObject("hZIn"));
119 sprintf(title,"Z coordinate exiting: %s ",mod->GetName());
120 list->Add(new TH1F("hZOut",strcat(title,mod->GetTitle()),100,
121 mod->ZMin(),mod->ZMax()));
122 dir->Remove(dir->FindObject("hZOut"));
125 gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
127 fDetDone = new Int_t[fNdets];
130 // Global QA histograms
132 fQAHist = new TObjArray(2);
133 fNvolumes=gMC->NofVolumes();
135 fQAHist->Add(new TH1F("hMCVcalls","Monte Carlo calls per volume",
136 fNvolumes, 0.5, fNvolumes+0.5));
137 h = (TH1F*) dir->FindObject("hMCVcalls");
138 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddVolumeName()"));
139 dir->Remove(dir->FindObject("hMCVcalls"));
141 // Build list of volume names
143 fVolNames=new TObjArray(fNvolumes);
144 for(i=0;i<fNvolumes;++i) {
145 AliModule *mod = (AliModule*)
146 (*gAlice->Modules())[gAlice->DetFromMate(gMC->VolId2Mate(i+1))];
147 (*fVolNames)[i]=new TNamed(gMC->VolName(i+1),mod->GetName());
150 fQAHist->Add(new TH1F("hMCMcalls","Monte Carlo calls per module",
151 fNdets, -0.5, fNdets-0.5));
152 h = (TH1F*) dir->FindObject("hMCMcalls");
153 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddModuleName()"));
155 dir->Remove(dir->FindObject("hMCMcalls"));
157 // Build list of module names
159 fModNames=new TObjArray(fNdets);
160 for(i=0;i<fNdets;++i)
162 new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
165 //_____________________________________________________________________________
167 AliMCQA::~AliMCQA() {
168 gROOT->GetListOfBrowsables()->Remove(this);
194 //_____________________________________________________________________________
195 void AliMCQA::Browse(TBrowser *b)
198 // Called when the item "Run" is clicked on the left pane
199 // of the Root browser.
200 // It displays the Root Trees and all detectors.
204 // Global histos first
206 TIter global(fQAHist);
207 while((hist = (TH1*)global()))
208 b->Add(hist,hist->GetTitle());
210 // Module histograms now
214 while((histos = (TList*)next())) {
216 while((hist = (TH1*)next1()))
217 b->Add(hist,hist->GetTitle());
221 //_____________________________________________________________________________
222 void AliMCQA::PreTrack()
225 // Called before each track
228 for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
231 //_____________________________________________________________________________
232 void AliMCQA::StepManager(Int_t id)
235 // Called at each step
240 // Fill Global histograms first
244 static TH1F* mcvcalls = (TH1F*) fQAHist->FindObject("hMCVcalls");
245 mcvcalls->Fill(gMC->CurrentVolID(copy));
246 static TH1F* mcmcalls = (TH1F*) fQAHist->FindObject("hMCMcalls");
250 // Now the step manager histograms
253 static Double_t mpi0=0;
254 static Double_t mpip=0;
255 static Double_t mpim=0;
256 static Double_t mep=0;
257 static Double_t mem=0;
259 Int_t num = gMC->TrackPid();
263 if (mpi0==0) mpi0=gAlice->PDGDB()->GetParticle(num)->Mass();
267 if (mpip==0) mpip=gAlice->PDGDB()->GetParticle(num)->Mass();
271 if (mpim==0) mpim=gAlice->PDGDB()->GetParticle(num)->Mass();
275 if (mep==0) mep=gAlice->PDGDB()->GetParticle(num)->Mass();
279 if (mem==0) mem=gAlice->PDGDB()->GetParticle(num)->Mass();
283 mass =gAlice->PDGDB()->GetParticle(num)->Mass();
287 static TLorentzVector p, x;
288 gMC->TrackMomentum(p);
289 gMC->TrackPosition(x);
290 Double_t energy = TMath::Max(p[3]-mass,1.e-12);
292 if(!fDetDone[fOldId] && !gMC->IsNewTrack()) {
293 TList *histold = (TList*) (*fQAList)[fOldId];
294 hist = (TH1F*) histold->FindObject("hEnOut");
295 hist->Fill(TMath::Log10(energy));
296 hist = (TH1F*) histold->FindObject("hZOut");
301 if(!fDetDone[id] && !gMC->IsNewTrack()) {
302 TList *histnew = (TList*) (*fQAList)[id];
303 hist = (TH1F*) histnew->FindObject("hEnIn");
304 hist->Fill(TMath::Log10(energy));
305 hist = (TH1F*) histnew->FindObject("hZIn");
312 //_____________________________________________________________________________
313 void AliMCQA::AddModuleName()
316 // Add function DrawModuleName to the module frequency histogram pad
318 static TVirtualPad* oldpad=0;
320 gPad->GetCanvas()->FeedbackMode(kTRUE);
321 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
322 DrawPaveLabel(fMPaveLabel);
327 //_____________________________________________________________________________
328 void AliMCQA::AddVolumeName()
331 // Add function DrawVolumeName to the volume frequency histogram pad
333 static TVirtualPad* oldpad=0;
335 gPad->GetCanvas()->FeedbackMode(kTRUE);
336 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
337 DrawPaveLabel(fVPaveLabel);
342 //_____________________________________________________________________________
343 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
346 // Draws the PaveLabel with the meaning of the bin
348 float uxmin = gPad->GetUxmin();
349 float uxmax = gPad->GetUxmax();
350 float uymin = gPad->GetUymin();
351 float uymax = gPad->GetUymax();
352 float lx = uxmax-uxmin;
353 float ly = uymax-uymin;
357 = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
358 uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
363 //_____________________________________________________________________________
364 Int_t AliMCQA::GetHBin(const char* hname)
367 // Get the bin where the cursor is
369 TList *dir = gDirectory->GetList();
370 TH1 *h=(TH1*)dir->FindObject(hname);
373 int px = gPad->GetEventX();
374 Float_t upx = gPad->AbsPixeltoX(px);
375 Float_t x = gPad->PadtoX(upx);
377 return h->GetXaxis()->FindBin(x);
380 //_____________________________________________________________________________
381 void AliMCQA::DrawModuleName()
384 // Writes the name of the module of the bin where the cursor is
386 TObject *select = gPad->GetSelected();
389 Int_t binx = GetHBin("hMCMcalls");
390 if(0<binx && binx<=fNdets) {
392 strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
393 fMPaveLabel->SetLabel(lab);
400 //_____________________________________________________________________________
401 void AliMCQA::DrawVolumeName()
404 // Writes the name of the volume:module of the bin where the cursor is
406 TObject *select = gPad->GetSelected();
409 Int_t binx = GetHBin("hMCVcalls");
410 if(0<binx && binx<=fNvolumes) {
412 sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
413 ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
414 fVPaveLabel->SetLabel(lab);