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.4 2001/01/17 10:50:50 hristov
19 Corrections to destructors
21 Revision 1.3 2000/12/18 14:16:31 alibrary
24 Revision 1.2 2000/12/18 11:33:48 alibrary
25 New call frequence histograms per module and volume
27 Revision 1.1 2000/11/30 07:12:48 alibrary
28 Introducing new Rndm and QA classes
32 ///////////////////////////////////////////////////////////////////////////////
35 ///////////////////////////////////////////////////////////////////////////////
40 #include "TObjArray.h"
46 #include "TLorentzVector.h"
47 #include "TDatabasePDG.h"
51 #include "TPaveLabel.h"
56 #include "AliModule.h"
62 //_____________________________________________________________________________
63 AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
64 fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
67 // Default constructor
71 //_____________________________________________________________________________
72 AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
75 // Constructor, creates the list of lists of histograms
83 fQAList = new TObjArray(ndets);
84 TObjArray &hist = *fQAList;
88 TObjArray &mods = *(gAlice->Modules());
90 TList *dir = gDirectory->GetList();
91 for (i=0; i<ndets; i++) {
92 hist[i] = list = new TList();
93 mod = (AliModule *) mods[i];
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"));
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"));
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"));
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"));
116 gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
118 fDetDone = new Int_t[fNdets];
121 // Global QA histograms
123 fQAHist = new TObjArray(2);
124 fNvolumes=gMC->NofVolumes();
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"));
132 // Build list of volume names
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());
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()"));
146 dir->Remove(dir->FindObject("hMCMcalls"));
148 // Build list of module names
150 fModNames=new TObjArray(fNdets);
151 for(i=0;i<fNdets;++i)
153 new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
156 //_____________________________________________________________________________
158 AliMCQA::~AliMCQA() {
181 //_____________________________________________________________________________
182 void AliMCQA::Browse(TBrowser *b)
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.
191 // Global histos first
193 TIter global(fQAHist);
194 while((hist = (TH1*)global()))
195 b->Add(hist,hist->GetTitle());
197 // Module histograms now
201 while((histos = (TList*)next())) {
203 while((hist = (TH1*)next1()))
204 b->Add(hist,hist->GetTitle());
208 //_____________________________________________________________________________
209 void AliMCQA::PreTrack()
212 // Called before each track
215 for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
218 //_____________________________________________________________________________
219 void AliMCQA::StepManager(Int_t id)
222 // Called at each step
227 // Fill Global histograms first
231 static TH1F* mcvcalls = (TH1F*) fQAHist->FindObject("hMCVcalls");
232 mcvcalls->Fill(gMC->CurrentVolID(copy));
233 static TH1F* mcmcalls = (TH1F*) fQAHist->FindObject("hMCMcalls");
237 // Now the step manager histograms
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;
246 Int_t num = gMC->TrackPid();
250 if (mpi0==0) mpi0=gAlice->PDGDB()->GetParticle(num)->Mass();
254 if (mpip==0) mpip=gAlice->PDGDB()->GetParticle(num)->Mass();
258 if (mpim==0) mpim=gAlice->PDGDB()->GetParticle(num)->Mass();
262 if (mep==0) mep=gAlice->PDGDB()->GetParticle(num)->Mass();
266 if (mem==0) mem=gAlice->PDGDB()->GetParticle(num)->Mass();
270 mass =gAlice->PDGDB()->GetParticle(num)->Mass();
274 static TLorentzVector p, x;
275 gMC->TrackMomentum(p);
276 gMC->TrackPosition(x);
277 Double_t energy = TMath::Max(p[3]-mass,1.e-12);
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");
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");
299 //_____________________________________________________________________________
300 void AliMCQA::AddModuleName()
303 // Add function DrawModuleName to the module frequency histogram pad
305 static TVirtualPad* oldpad=0;
307 gPad->GetCanvas()->FeedbackMode(kTRUE);
308 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
309 DrawPaveLabel(fMPaveLabel);
314 //_____________________________________________________________________________
315 void AliMCQA::AddVolumeName()
318 // Add function DrawVolumeName to the volume frequency histogram pad
320 static TVirtualPad* oldpad=0;
322 gPad->GetCanvas()->FeedbackMode(kTRUE);
323 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
324 DrawPaveLabel(fVPaveLabel);
329 //_____________________________________________________________________________
330 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
333 // Draws the PaveLabel with the meaning of the bin
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;
344 = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
345 uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
350 //_____________________________________________________________________________
351 Int_t AliMCQA::GetHBin(const char* hname)
354 // Get the bin where the cursor is
356 TList *dir = gDirectory->GetList();
357 TH1 *h=(TH1*)dir->FindObject(hname);
360 int px = gPad->GetEventX();
361 Float_t upx = gPad->AbsPixeltoX(px);
362 Float_t x = gPad->PadtoX(upx);
364 return h->GetXaxis()->FindBin(x);
367 //_____________________________________________________________________________
368 void AliMCQA::DrawModuleName()
371 // Writes the name of the module of the bin where the cursor is
373 TObject *select = gPad->GetSelected();
376 Int_t binx = GetHBin("hMCMcalls");
377 if(0<binx && binx<=fNdets) {
379 strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
380 fMPaveLabel->SetLabel(lab);
387 //_____________________________________________________________________________
388 void AliMCQA::DrawVolumeName()
391 // Writes the name of the volume:module of the bin where the cursor is
393 TObject *select = gPad->GetSelected();
396 Int_t binx = GetHBin("hMCVcalls");
397 if(0<binx && binx<=fNvolumes) {
399 sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
400 ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
401 fVPaveLabel->SetLabel(lab);