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.3 2000/12/18 14:16:31 alibrary
21 Revision 1.2 2000/12/18 11:33:48 alibrary
22 New call frequence histograms per module and volume
24 Revision 1.1 2000/11/30 07:12:48 alibrary
25 Introducing new Rndm and QA classes
29 ///////////////////////////////////////////////////////////////////////////////
32 ///////////////////////////////////////////////////////////////////////////////
37 #include "TObjArray.h"
43 #include "TLorentzVector.h"
44 #include "TDatabasePDG.h"
48 #include "TPaveLabel.h"
53 #include "AliModule.h"
59 //_____________________________________________________________________________
60 AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
61 fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
64 // Default constructor
68 //_____________________________________________________________________________
69 AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
72 // Constructor, creates the list of lists of histograms
80 fQAList = new TObjArray(ndets);
81 TObjArray &hist = *fQAList;
85 TObjArray &mods = *(gAlice->Modules());
87 TList *dir = gDirectory->GetList();
88 for (i=0; i<ndets; i++) {
89 hist[i] = list = new TList();
90 mod = (AliModule *) mods[i];
93 sprintf(title,"Spectrum entering: %s ",mod->GetName());
94 list->Add(new TH1F("hEnIn",strcat(title,mod->GetTitle()),100,-4,2));
95 dir->Remove(dir->FindObject("hEnIn"));
97 sprintf(title,"Spectrum exiting: %s ",mod->GetName());
98 list->Add(new TH1F("hEnOut",strcat(title,mod->GetTitle()),100,-4,2));
99 dir->Remove(dir->FindObject("hEnOut"));
102 sprintf(title,"Z coordinate entering: %s ",mod->GetName());
103 list->Add(new TH1F("hZIn",strcat(title,mod->GetTitle()),100,
104 mod->ZMin(),mod->ZMax()));
105 dir->Remove(dir->FindObject("hZIn"));
107 sprintf(title,"Z coordinate exiting: %s ",mod->GetName());
108 list->Add(new TH1F("hZOut",strcat(title,mod->GetTitle()),100,
109 mod->ZMin(),mod->ZMax()));
110 dir->Remove(dir->FindObject("hZOut"));
113 gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
115 fDetDone = new Int_t[fNdets];
118 // Global QA histograms
120 fQAHist = new TObjArray(2);
121 fNvolumes=gMC->NofVolumes();
123 fQAHist->Add(new TH1F("hMCVcalls","Monte Carlo calls per volume",
124 fNvolumes, 0.5, fNvolumes+0.5));
125 h = (TH1F*) dir->FindObject("hMCVcalls");
126 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddVolumeName()"));
127 dir->Remove(dir->FindObject("hMCVcalls"));
129 // Build list of volume names
131 fVolNames=new TObjArray(fNvolumes);
132 for(i=0;i<fNvolumes;++i) {
133 AliModule *mod = (AliModule*)
134 (*gAlice->Modules())[gAlice->DetFromMate(gMC->VolId2Mate(i+1))];
135 (*fVolNames)[i]=new TNamed(gMC->VolName(i+1),mod->GetName());
138 fQAHist->Add(new TH1F("hMCMcalls","Monte Carlo calls per module",
139 fNdets, -0.5, fNdets-0.5));
140 h = (TH1F*) dir->FindObject("hMCMcalls");
141 h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddModuleName()"));
143 dir->Remove(dir->FindObject("hMCMcalls"));
145 // Build list of module names
147 fModNames=new TObjArray(fNdets);
148 for(i=0;i<fNdets;++i)
150 new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
153 //_____________________________________________________________________________
155 AliMCQA::~AliMCQA() {
178 //_____________________________________________________________________________
179 void AliMCQA::Browse(TBrowser *b)
182 // Called when the item "Run" is clicked on the left pane
183 // of the Root browser.
184 // It displays the Root Trees and all detectors.
188 // Global histos first
190 TIter global(fQAHist);
191 while((hist = (TH1*)global()))
192 b->Add(hist,hist->GetTitle());
194 // Module histograms now
198 while((histos = (TList*)next())) {
200 while((hist = (TH1*)next1()))
201 b->Add(hist,hist->GetTitle());
205 //_____________________________________________________________________________
206 void AliMCQA::PreTrack()
209 // Called before each track
212 for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
215 //_____________________________________________________________________________
216 void AliMCQA::StepManager(Int_t id)
219 // Called at each step
224 // Fill Global histograms first
226 hist = (TH1F*) fQAHist->FindObject("hMCVcalls");
227 hist->Fill(gMC->CurrentVolID(copy));
228 hist = (TH1F*) fQAHist->FindObject("hMCMcalls");
231 // Now the step manager histograms
235 gMC->TrackMomentum(p);
236 gMC->TrackPosition(x);
237 Double_t energy = TMath::Max(
238 p[3]-gAlice->PDGDB()->GetParticle(gMC->TrackPid())->Mass(),1.e-12);
240 if(!fDetDone[fOldId] && !gMC->IsNewTrack()) {
241 TList *histold = (TList*) (*fQAList)[fOldId];
242 hist = (TH1F*) histold->FindObject("hEnOut");
243 hist->Fill(TMath::Log10(energy));
244 hist = (TH1F*) histold->FindObject("hZOut");
249 if(!fDetDone[id] && !gMC->IsNewTrack()) {
250 TList *histnew = (TList*) (*fQAList)[id];
251 hist = (TH1F*) histnew->FindObject("hEnIn");
252 hist->Fill(TMath::Log10(energy));
253 hist = (TH1F*) histnew->FindObject("hZIn");
260 //_____________________________________________________________________________
261 void AliMCQA::AddModuleName()
264 // Add function DrawModuleName to the module frequency histogram pad
266 static TVirtualPad* oldpad=0;
268 gPad->GetCanvas()->FeedbackMode(kTRUE);
269 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
270 DrawPaveLabel(fMPaveLabel);
275 //_____________________________________________________________________________
276 void AliMCQA::AddVolumeName()
279 // Add function DrawVolumeName to the volume frequency histogram pad
281 static TVirtualPad* oldpad=0;
283 gPad->GetCanvas()->FeedbackMode(kTRUE);
284 if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
285 DrawPaveLabel(fVPaveLabel);
290 //_____________________________________________________________________________
291 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
294 // Draws the PaveLabel with the meaning of the bin
296 float uxmin = gPad->GetUxmin();
297 float uxmax = gPad->GetUxmax();
298 float uymin = gPad->GetUymin();
299 float uymax = gPad->GetUymax();
300 float lx = uxmax-uxmin;
301 float ly = uymax-uymin;
305 = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
306 uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
311 //_____________________________________________________________________________
312 Int_t AliMCQA::GetHBin(const char* hname)
315 // Get the bin where the cursor is
317 TList *dir = gDirectory->GetList();
318 TH1 *h=(TH1*)dir->FindObject(hname);
321 int px = gPad->GetEventX();
322 Float_t upx = gPad->AbsPixeltoX(px);
323 Float_t x = gPad->PadtoX(upx);
325 return h->GetXaxis()->FindBin(x);
328 //_____________________________________________________________________________
329 void AliMCQA::DrawModuleName()
332 // Writes the name of the module of the bin where the cursor is
334 TObject *select = gPad->GetSelected();
337 Int_t binx = GetHBin("hMCMcalls");
338 if(0<binx && binx<=fNdets) {
340 strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
341 fMPaveLabel->SetLabel(lab);
348 //_____________________________________________________________________________
349 void AliMCQA::DrawVolumeName()
352 // Writes the name of the volume:module of the bin where the cursor is
354 TObject *select = gPad->GetSelected();
357 Int_t binx = GetHBin("hMCVcalls");
358 if(0<binx && binx<=fNvolumes) {
360 sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
361 ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
362 fVPaveLabel->SetLabel(lab);