]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCQA.cxx
Merging the VirtualMC branch to the main development branch (HEAD)
[u/mrichter/AliRoot.git] / STEER / AliMCQA.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.6.6.1  2002/10/12 21:41:00  hristov
19 Remove the object from the list of browsables
20
21 Revision 1.6  2001/12/05 08:31:25  hristov
22 Destructor corrected, thanks to R.Brun
23
24 Revision 1.5  2001/01/26 19:58:48  hristov
25 Major upgrade of AliRoot code
26
27 Revision 1.4  2001/01/17 10:50:50  hristov
28 Corrections to destructors
29
30 Revision 1.3  2000/12/18 14:16:31  alibrary
31 HP compatibility fix
32
33 Revision 1.2  2000/12/18 11:33:48  alibrary
34 New call frequence histograms per module and volume
35
36 Revision 1.1  2000/11/30 07:12:48  alibrary
37 Introducing new Rndm and QA classes
38
39 */
40
41 ///////////////////////////////////////////////////////////////////////////////
42 //                                                                           //
43 //                                                                           //
44 ///////////////////////////////////////////////////////////////////////////////
45
46
47 #include <strings.h>
48
49 #include "TObjArray.h"
50 #include "TH1.h"
51 #include "TList.h"
52 #include "TROOT.h"
53 #include "TBrowser.h"
54 #include "TMath.h"
55 #include "TLorentzVector.h"
56 #include "TDatabasePDG.h"
57 #include "TMath.h"
58 #include "TPad.h"
59 #include "TExec.h"
60 #include "TPaveLabel.h"
61 #include "TCanvas.h"
62
63 #include "AliMCQA.h"
64 #include "AliRun.h"
65 #include "AliModule.h"
66 #include "AliMC.h"
67
68 ClassImp(AliMCQA)
69
70
71 //_____________________________________________________________________________
72 AliMCQA::AliMCQA() : fQAList(0), fDetDone(0), fQAHist(0), fVolNames(0),
73                      fModNames(0),fMPaveLabel(0),fVPaveLabel(0)
74 {
75   //
76   // Default constructor
77   //
78 }
79
80 //_____________________________________________________________________________
81 AliMCQA::AliMCQA(Int_t ndets) : fMPaveLabel(0),fVPaveLabel(0)
82 {
83   //
84   // Constructor, creates the list of lists of histograms
85   //
86   TList *list;
87   TH1F* h;
88   Int_t i;
89   
90   fNdets=ndets;
91
92   fQAList = new TObjArray(ndets);
93   TObjArray &hist = *fQAList;
94
95   char title[100];
96   //
97   TObjArray &mods = *(gAlice->Modules());
98   AliModule *mod;
99   TList *dir = gDirectory->GetList();
100   for (i=0; i<ndets; i++) {
101     hist[i] = list = new TList();
102     mod = (AliModule *) mods[i];
103
104     // Energy Spectrum
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"));
108
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"));
112
113     // Z position
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"));
118
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"));
123   }
124   //
125   gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
126
127   fDetDone = new Int_t[fNdets];
128
129   //
130   // Global QA histograms
131   //
132   fQAHist = new TObjArray(2);
133   fNvolumes=gMC->NofVolumes();
134   
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"));
140   //
141   // Build list of volume names
142   //
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());
148   }
149
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()"));
154
155   dir->Remove(dir->FindObject("hMCMcalls"));
156   //
157   // Build list of module names
158   //
159   fModNames=new TObjArray(fNdets);
160   for(i=0;i<fNdets;++i) 
161     (*fModNames)[i]=
162       new TNamed(((AliModule *)(*gAlice->Modules())[i])->GetName(),"");
163 }
164
165 //_____________________________________________________________________________
166
167 AliMCQA::~AliMCQA() {
168   gROOT->GetListOfBrowsables()->Remove(this);
169   if (fQAList) {
170     fQAList->Delete();
171     delete fQAList;
172     fQAList=0;
173   }
174   if (fQAHist) {
175     fQAHist->Delete();
176     delete fQAHist;
177     fQAHist=0;
178   }
179   if (fVolNames) {
180     fVolNames->Delete();
181     delete fVolNames;
182     fVolNames=0;
183   }
184   if (fModNames) {
185     fModNames->Delete();
186     delete fModNames;
187     fModNames=0;
188   }
189   delete [] fDetDone;
190   delete fMPaveLabel;
191   delete fVPaveLabel;
192 }
193
194 //_____________________________________________________________________________
195 void AliMCQA::Browse(TBrowser *b)
196 {
197   //
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.
201   //
202   TH1 *hist;
203   //
204   // Global histos first
205   //
206   TIter global(fQAHist);
207   while((hist = (TH1*)global())) 
208     b->Add(hist,hist->GetTitle());
209   //
210   // Module histograms now
211   //
212   TIter next(fQAList);
213   TList *histos;
214   while((histos = (TList*)next())) {
215     TIter next1(histos);
216     while((hist = (TH1*)next1())) 
217       b->Add(hist,hist->GetTitle());
218   }
219 }
220
221 //_____________________________________________________________________________
222 void AliMCQA::PreTrack()
223 {
224   //
225   // Called before each track
226   //
227   fOldId=-1;
228   for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
229 }
230
231 //_____________________________________________________________________________
232 void AliMCQA::StepManager(Int_t id)
233 {
234   //
235   // Called at each step
236   //
237   TH1F *hist;
238   Int_t copy;
239   //
240   // Fill Global histograms first
241   //
242
243
244   static TH1F* mcvcalls = (TH1F*) fQAHist->FindObject("hMCVcalls");
245   mcvcalls->Fill(gMC->CurrentVolID(copy));
246   static TH1F* mcmcalls = (TH1F*) fQAHist->FindObject("hMCMcalls");
247   mcmcalls->Fill(id);
248
249   //
250   // Now the step manager histograms
251   //
252   if(fOldId != id) {
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;
258     Double_t mass=0;
259     Int_t num = gMC->TrackPid();
260
261     switch (num) {
262     case 111: 
263       if (mpi0==0) mpi0=gAlice->PDGDB()->GetParticle(num)->Mass(); 
264       mass=mpi0;
265       break;
266     case 211:
267       if (mpip==0) mpip=gAlice->PDGDB()->GetParticle(num)->Mass(); 
268       mass=mpip;
269       break;
270     case -211:
271       if (mpim==0) mpim=gAlice->PDGDB()->GetParticle(num)->Mass(); 
272       mass=mpim;
273       break;
274     case 11:
275       if (mep==0) mep=gAlice->PDGDB()->GetParticle(num)->Mass(); 
276       mass=mep;
277       break;
278     case -11:
279       if (mem==0) mem=gAlice->PDGDB()->GetParticle(num)->Mass(); 
280       mass=mem;
281       break;
282     default:
283       mass =gAlice->PDGDB()->GetParticle(num)->Mass();
284       break; 
285     }
286
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);
291     if(fOldId > -1) {
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");
297         hist->Fill(x[2]);
298         fDetDone[fOldId]=1;
299       }
300     }
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");
306       hist->Fill(x[2]);
307     }
308     fOldId=id;
309   }
310 }
311
312 //_____________________________________________________________________________
313 void AliMCQA::AddModuleName()
314 {
315   //
316   // Add function DrawModuleName to the module frequency histogram pad
317   //
318   static TVirtualPad* oldpad=0;
319   if(oldpad!=gPad) {
320     gPad->GetCanvas()->FeedbackMode(kTRUE);
321     if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
322     DrawPaveLabel(fMPaveLabel);
323     oldpad=gPad;
324   }
325 }
326
327 //_____________________________________________________________________________
328 void AliMCQA::AddVolumeName()
329 {
330   //
331   // Add function DrawVolumeName to the volume frequency histogram pad
332   //
333   static TVirtualPad* oldpad=0;
334   if(oldpad!=gPad) {
335     gPad->GetCanvas()->FeedbackMode(kTRUE);
336     if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
337     DrawPaveLabel(fVPaveLabel);
338     oldpad=gPad;
339   }
340 }
341
342 //_____________________________________________________________________________
343 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
344 {
345   //
346   // Draws the PaveLabel with the meaning of the bin
347   //
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;
354
355   if(pv) delete pv;
356   pv 
357     = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
358                      uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
359                      "");
360   pv->Draw();
361 }
362
363 //_____________________________________________________________________________
364 Int_t AliMCQA::GetHBin(const char* hname)
365 {
366   //
367   // Get the bin where the cursor is
368   //
369   TList *dir = gDirectory->GetList();
370   TH1 *h=(TH1*)dir->FindObject(hname);
371   
372
373   int px = gPad->GetEventX();
374   Float_t upx = gPad->AbsPixeltoX(px);
375   Float_t x = gPad->PadtoX(upx);
376     
377   return h->GetXaxis()->FindBin(x);
378 }
379
380 //_____________________________________________________________________________
381 void AliMCQA::DrawModuleName()
382 {
383   //
384   // Writes the name of the module of the bin where the cursor is
385   //
386   TObject *select = gPad->GetSelected();
387   if(!select) return;
388   
389   Int_t binx = GetHBin("hMCMcalls");
390   if(0<binx && binx<=fNdets) {
391     char lab[15];
392     strcpy(lab,((TNamed*)(*fModNames)[binx-1])->GetName());
393     fMPaveLabel->SetLabel(lab);
394   
395     gPad->Modified();
396     gPad->Update();
397   }
398 }
399
400 //_____________________________________________________________________________
401 void AliMCQA::DrawVolumeName()
402 {
403   //
404   // Writes the name of the volume:module of the bin where the cursor is
405   //
406   TObject *select = gPad->GetSelected();
407   if(!select) return;
408
409   Int_t binx = GetHBin("hMCVcalls");
410   if(0<binx && binx<=fNvolumes) {
411     char lab[20];
412     sprintf(lab,"%s: %s",((TNamed*)(*fVolNames)[binx-1])->GetName(),
413             ((TNamed*)(*fVolNames)[binx-1])->GetTitle());
414     fVPaveLabel->SetLabel(lab);
415     
416     gPad->Modified();
417     gPad->Update();
418   }
419 }
420
421