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.2 2000/12/18 11:33:48 alibrary
19 New call frequence histograms per module and volume
21 Revision 1.1 2000/11/30 07:12:48 alibrary
22 Introducing new Rndm and QA classes
26 ///////////////////////////////////////////////////////////////////////////////
29 ///////////////////////////////////////////////////////////////////////////////
34 #include "TObjArray.h"
40 #include "TLorentzVector.h"
41 #include "TDatabasePDG.h"
45 #include "TPaveLabel.h"
50 #include "AliModule.h"
56 //_____________________________________________________________________________
57 AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
58 fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
61 // Default constructor
65 //_____________________________________________________________________________
66 AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
69 // Constructor, creates the list of lists of histograms
77 fQAList = new TObjArray(ndets);
78 TObjArray &hist = *fQAList;
82 TObjArray &mods = *(gAlice->Modules());
84 TList *dir = gDirectory->GetList();
85 for (i=0; i<ndets; i++) {
86 hist[i] = list = new TList();
87 mod = (AliModule *) mods[i];
90 sprintf(title,"Spectrum entering: %s ",mod->GetName());
91 list->Add(new TH1F("hEnIn",strcat(title,mod->GetTitle()),100,-4,2));
92 dir->Remove(dir->FindObject("hEnIn"));
94 sprintf(title,"Spectrum exiting: %s ",mod->GetName());
95 list->Add(new TH1F("hEnOut",strcat(title,mod->GetTitle()),100,-4,2));
96 dir->Remove(dir->FindObject("hEnOut"));
99 sprintf(title,"Z coordinate entering: %s ",mod->GetName());
100 list->Add(new TH1F("hZIn",strcat(title,mod->GetTitle()),100,
101 mod->ZMin(),mod->ZMax()));
102 dir->Remove(dir->FindObject("hZIn"));
104 sprintf(title,"Z coordinate exiting: %s ",mod->GetName());
105 list->Add(new TH1F("hZOut",strcat(title,mod->GetTitle()),100,
106 mod->ZMin(),mod->ZMax()));
107 dir->Remove(dir->FindObject("hZOut"));
110 gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
112 fDetDone = new Int_t[fNdets];
115 // Global QA histograms
117 fQAHist = new TObjArray(2);
118 fNvolumes=gMC->NofVolumes();
120 fQAHist->Add(new TH1F("hMCVcalls","Monte Carlo calls per volume",
121 fNvolumes, 0.5, fNvolumes+0.5));
122 h = (TH1F*) dir->FindObject("hMCVcalls");
123 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddVolumeName()"));
124 dir->Remove(dir->FindObject("hMCVcalls"));
126 // Build list of volume names
128 fVolNames=new TObjArray(fNvolumes);
129 for(i=0;i<fNvolumes;++i) {
130 AliModule *mod = (AliModule*)
131 (*gAlice->Modules())[gAlice->DetFromMate(gMC->VolId2Mate(i+1))];
132 (*fVolNames)[i]=new TNamed(gMC->VolName(i+1),mod->GetName());
135 fQAHist->Add(new TH1F("hMCMcalls","Monte Carlo calls per module",
136 fNdets, -0.5, fNdets-0.5));
137 h = (TH1F*) dir->FindObject("hMCMcalls");
138 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddModuleName()"));
140 dir->Remove(dir->FindObject("hMCMcalls"));
142 // Build list of module names
144 fModNames=new TObjArray(fNdets);
145 for(i=0;i<fNdets;++i)
147 new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
150 //_____________________________________________________________________________
151 void AliMCQA::Browse(TBrowser *b)
154 // Called when the item "Run" is clicked on the left pane
155 // of the Root browser.
156 // It displays the Root Trees and all detectors.
160 // Global histos first
162 TIter global(fQAHist);
163 while((hist = (TH1*)global()))
164 b->Add(hist,hist->GetTitle());
166 // Module histograms now
170 while((histos = (TList*)next())) {
172 while((hist = (TH1*)next1()))
173 b->Add(hist,hist->GetTitle());
177 //_____________________________________________________________________________
178 void AliMCQA::PreTrack()
181 // Called before each track
184 for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
187 //_____________________________________________________________________________
188 void AliMCQA::StepManager(Int_t id)
191 // Called at each step
196 // Fill Global histograms first
198 hist = (TH1F*) fQAHist->FindObject("hMCVcalls");
199 hist->Fill(gMC->CurrentVolID(copy));
200 hist = (TH1F*) fQAHist->FindObject("hMCMcalls");
203 // Now the step manager histograms
207 gMC->TrackMomentum(p);
208 gMC->TrackPosition(x);
209 Double_t energy = TMath::Max(
210 p[3]-gAlice->PDGDB()->GetParticle(gMC->TrackPid())->Mass(),1.e-12);
212 if(!fDetDone[fOldId] && !gMC->IsNewTrack()) {
213 TList *histold = (TList*) (*fQAList)[fOldId];
214 hist = (TH1F*) histold->FindObject("hEnOut");
215 hist->Fill(TMath::Log10(energy));
216 hist = (TH1F*) histold->FindObject("hZOut");
221 if(!fDetDone[id] && !gMC->IsNewTrack()) {
222 TList *histnew = (TList*) (*fQAList)[id];
223 hist = (TH1F*) histnew->FindObject("hEnIn");
224 hist->Fill(TMath::Log10(energy));
225 hist = (TH1F*) histnew->FindObject("hZIn");
232 //_____________________________________________________________________________
233 void AliMCQA::AddModuleName()
236 // Add function DrawModuleName to the module frequency histogram pad
238 static TVirtualPad* oldpad=0;
240 gPad->GetCanvas()->FeedbackMode(kTRUE);
241 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
242 DrawPaveLabel(fMPaveLabel);
247 //_____________________________________________________________________________
248 void AliMCQA::AddVolumeName()
251 // Add function DrawVolumeName to the volume frequency histogram pad
253 static TVirtualPad* oldpad=0;
255 gPad->GetCanvas()->FeedbackMode(kTRUE);
256 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
257 DrawPaveLabel(fVPaveLabel);
262 //_____________________________________________________________________________
263 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
266 // Draws the PaveLabel with the meaning of the bin
268 float uxmin = gPad->GetUxmin();
269 float uxmax = gPad->GetUxmax();
270 float uymin = gPad->GetUymin();
271 float uymax = gPad->GetUymax();
272 float lx = uxmax-uxmin;
273 float ly = uymax-uymin;
277 = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
278 uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
283 //_____________________________________________________________________________
284 Int_t AliMCQA::GetHBin(const char* hname)
287 // Get the bin where the cursor is
289 TList *dir = gDirectory->GetList();
290 TH1 *h=(TH1*)dir->FindObject(hname);
293 int px = gPad->GetEventX();
294 Float_t upx = gPad->AbsPixeltoX(px);
295 Float_t x = gPad->PadtoX(upx);
297 return h->GetXaxis()->FindBin(x);
300 //_____________________________________________________________________________
301 void AliMCQA::DrawModuleName()
304 // Writes the name of the module of the bin where the cursor is
306 TObject *select = gPad->GetSelected();
309 Int_t binx = GetHBin("hMCMcalls");
310 if(0<binx && binx<=fNdets) {
312 strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
313 fMPaveLabel->SetLabel(lab);
320 //_____________________________________________________________________________
321 void AliMCQA::DrawVolumeName()
324 // Writes the name of the volume:module of the bin where the cursor is
326 TObject *select = gPad->GetSelected();
329 Int_t binx = GetHBin("hMCVcalls");
330 if(0<binx && binx<=fNvolumes) {
332 sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
333 ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
334 fVPaveLabel->SetLabel(lab);