]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/AliMCQA.cxx
Coding violations
[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 /* $Id$ */
17
18 ///////////////////////////////////////////////////////////////////////////////
19 //                                                                           //
20 // Quality Assurance class                                                   //
21 // Provides standard histograms for every element of the detector            //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25
26 #include <strings.h>
27
28 #include <TBrowser.h>
29 #include <TCanvas.h>
30 #include <TDatabasePDG.h>
31 #include <TExec.h>
32 #include <TH1.h>
33 #include <TList.h>
34 #include <TLorentzVector.h>
35 #include <TMath.h>
36 #include <TObjArray.h>
37 #include <TPaveLabel.h>
38 #include <TROOT.h>
39 #include <TVirtualMC.h>
40
41 #include "AliMC.h"
42 #include "AliMCQA.h"
43 #include "AliModule.h"
44 #include "AliRun.h"
45
46 ClassImp(AliMCQA)
47
48 //_______________________________________________________________________
49 AliMCQA::AliMCQA():
50   fNdets(0),
51   fNvolumes(0),
52   fQAList(0),
53   fOldId(0),
54   fDetDone(0),
55   fQAHist(0),
56   fVolNames(0),
57   fModNames(0),
58   fMPaveLabel(0),
59   fVPaveLabel(0)
60 {
61   //
62   // Default constructor
63   //
64 }
65
66 //_______________________________________________________________________
67 AliMCQA::AliMCQA(const AliMCQA &qa):
68   TObject(qa),
69   fNdets(0),
70   fNvolumes(0),
71   fQAList(0),
72   fOldId(0),
73   fDetDone(0),
74   fQAHist(0),
75   fVolNames(0),
76   fModNames(0),
77   fMPaveLabel(0),
78   fVPaveLabel(0)
79 {
80   //
81   // Copy constructor
82   //
83   qa.Copy(*this);
84 }
85
86 //_______________________________________________________________________
87 AliMCQA::AliMCQA(Int_t ndets):
88   fNdets(ndets),
89   fNvolumes(gMC->NofVolumes()),
90   fQAList(new TObjArray(ndets)),
91   fOldId(0),
92   fDetDone(new Int_t[ndets]),
93   fQAHist(new TObjArray(2)),
94   fVolNames(new TObjArray(fNvolumes)),
95   fModNames(new TObjArray(fNdets)),
96   fMPaveLabel(0),
97   fVPaveLabel(0)
98 {
99   //
100   // Constructor, creates the list of lists of histograms
101   //
102   TList *list;
103   TH1F* h;
104   Int_t i;
105   
106   TObjArray &hist = *fQAList;
107
108   char title[100];
109   //
110   TObjArray &mods = *(gAlice->Modules());
111   TList *dir = gDirectory->GetList();
112   for (i=0; i<ndets; i++) {
113     hist[i] = list = new TList();
114     AliModule *mod = dynamic_cast<AliModule*>(mods[i]);
115
116     // Energy Spectrum
117     sprintf(title,"Spectrum entering: %s ",mod->GetName());
118     list->Add(new TH1F("hEnIn",strcat(title,mod->GetTitle()),100,-4,2));
119     dir->Remove(dir->FindObject("hEnIn"));
120
121     sprintf(title,"Spectrum exiting: %s ",mod->GetName());
122     list->Add(new TH1F("hEnOut",strcat(title,mod->GetTitle()),100,-4,2));
123     dir->Remove(dir->FindObject("hEnOut"));
124
125     // Z position
126     sprintf(title,"Z coordinate entering: %s ",mod->GetName());
127     list->Add(new TH1F("hZIn",strcat(title,mod->GetTitle()),100,
128                        mod->ZMin(),mod->ZMax()));
129     dir->Remove(dir->FindObject("hZIn"));
130
131     sprintf(title,"Z coordinate exiting: %s ",mod->GetName());
132     list->Add(new TH1F("hZOut",strcat(title,mod->GetTitle()),100,
133                        mod->ZMin(),mod->ZMax()));
134     dir->Remove(dir->FindObject("hZOut"));
135   }
136   //
137   gROOT->GetListOfBrowsables()->Add(this,"AliMCQA");
138
139   //
140   // Global QA histograms
141   //
142   
143   fQAHist->Add(new TH1F("hMCVcalls","Monte Carlo calls per volume",
144                         fNvolumes, 0.5, fNvolumes+0.5));
145   h = dynamic_cast<TH1F*>(dir->FindObject("hMCVcalls"));
146   h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddVolumeName()"));
147   dir->Remove(dir->FindObject("hMCVcalls"));
148   //
149   // Build list of volume names
150   //
151   AliMC * mc = gAlice->GetMCApp();
152   for(i=0;i<fNvolumes;++i) {
153     AliModule *mod = dynamic_cast<AliModule*>
154       ((*gAlice->Modules())[mc->DetFromMate(gMC->VolId2Mate(i+1))]);
155     (*fVolNames)[i]=new TNamed(gMC->VolName(i+1),mod->GetName());
156   }
157
158   fQAHist->Add(new TH1F("hMCMcalls","Monte Carlo calls per module",
159                         fNdets, -0.5, fNdets-0.5));
160   h = dynamic_cast<TH1F*>(dir->FindObject("hMCMcalls"));
161    h->GetListOfFunctions()->Add(new TExec("ex","gAlice->GetMCQA()->AddModuleName()"));
162
163   dir->Remove(dir->FindObject("hMCMcalls"));
164   //
165   // Build list of module names
166   //
167   for(i=0;i<fNdets;++i) 
168     (*fModNames)[i]=
169       new TNamed((dynamic_cast<AliModule*>((*gAlice->Modules())[i]))->GetName(),"");
170 }
171
172 //_______________________________________________________________________
173 void AliMCQA::Copy(TObject &) const
174 {
175   Fatal("Copy ctor","Not implemented!\n");
176 }
177
178 //_______________________________________________________________________
179 AliMCQA::~AliMCQA() 
180  {
181   //
182   // Destructor
183   //
184   gROOT->GetListOfBrowsables()->Remove(this);
185   //if program crashes here - it probobly means that 
186   //one of added browsables was deleted and not removed previously from that list
187   //skowron
188   
189   if (fQAList) {
190     fQAList->Delete();
191     delete fQAList;
192     fQAList = 0;
193   }
194   if (fQAHist) {
195     fQAHist->Delete();
196     delete fQAHist;
197     fQAHist = 0;
198   }
199   if (fVolNames) {
200     fVolNames->Delete();
201     delete fVolNames;
202     fVolNames = 0;
203   }
204   if (fModNames) {
205     fModNames->Delete();
206     delete fModNames;
207     fModNames = 0;
208   }
209   delete [] fDetDone;
210   delete fMPaveLabel;
211   delete fVPaveLabel;
212 }
213
214 //_______________________________________________________________________
215 void AliMCQA::Browse(TBrowser *b)
216 {
217   //
218   // Called when the item "Run" is clicked on the left pane
219   // of the Root browser.
220   // It displays the Root Trees and all detectors.
221   //
222   TH1 *hist;
223   //
224   // Global histos first
225   //
226   TIter global(fQAHist);
227   while((hist = dynamic_cast<TH1*>(global())))
228     b->Add(hist,hist->GetTitle());
229   //
230   // Module histograms now
231   //
232   TIter next(fQAList);
233   TList *histos;
234   while((histos = dynamic_cast<TList*>(next()))) {
235     TIter next1(histos);
236     while((hist = dynamic_cast<TH1*>(next1())))
237       b->Add(hist,hist->GetTitle());
238   }
239 }
240
241 //_______________________________________________________________________
242 void AliMCQA::PreTrack()
243 {
244   //
245   // Called before each track
246   //
247   fOldId=-1;
248   for(Int_t i=0; i<fNdets; i++) fDetDone[i]=0;
249 }
250
251 //_______________________________________________________________________
252 void AliMCQA::StepManager(Int_t id)
253 {
254   //
255   // Called at each step
256   //
257   TH1F *hist;
258   Int_t copy;
259   //
260   // Fill Global histograms first
261   //
262
263
264   static TH1F* mcvcalls = dynamic_cast<TH1F*>(fQAHist->FindObject("hMCVcalls"));
265   mcvcalls->Fill(gMC->CurrentVolID(copy));
266   static TH1F* mcmcalls = dynamic_cast<TH1F*>(fQAHist->FindObject("hMCMcalls"));
267   mcmcalls->Fill(id);
268
269   //
270   // Now the step manager histograms
271   //
272   if(fOldId != id) {
273     static  Double_t mpi0=0;
274     static  Double_t mpip=0;
275     static  Double_t mpim=0;
276     static  Double_t mep=0;
277     static  Double_t mem=0;
278     Double_t mass = 0;
279     Int_t num = gMC->TrackPid();
280
281     switch (num) {
282     case 111: 
283       if (mpi0==0) mpi0=gAlice->PDGDB()->GetParticle(num)->Mass(); 
284       mass=mpi0;
285       break;
286     case 211:
287       if (mpip==0) mpip=gAlice->PDGDB()->GetParticle(num)->Mass(); 
288       mass=mpip;
289       break;
290     case -211:
291       if (mpim==0) mpim=gAlice->PDGDB()->GetParticle(num)->Mass(); 
292       mass=mpim;
293       break;
294     case 11:
295       if (mep==0) mep=gAlice->PDGDB()->GetParticle(num)->Mass(); 
296       mass=mep;
297       break;
298     case -11:
299       if (mem==0) mem=gAlice->PDGDB()->GetParticle(num)->Mass(); 
300       mass=mem;
301       break;
302     default:
303         if (gAlice->PDGDB()->GetParticle(num))
304             mass = gAlice->PDGDB()->GetParticle(num)->Mass();
305       break; 
306     }
307
308     static TLorentzVector p, x;
309     gMC->TrackMomentum(p);
310     gMC->TrackPosition(x);
311     Double_t energy = TMath::Max(p[3]-mass,1.e-12);
312     if(fOldId > -1) {
313       if(!fDetDone[fOldId] && !gMC->IsNewTrack()) {
314         TList *histold = dynamic_cast<TList*>((*fQAList)[fOldId]);
315         hist = dynamic_cast<TH1F*>(histold->FindObject("hEnOut"));
316         hist->Fill(TMath::Log10(energy));
317         hist = dynamic_cast<TH1F*>(histold->FindObject("hZOut"));
318         hist->Fill(x[2]);
319         fDetDone[fOldId]=1;
320       }
321     }
322     if(!fDetDone[id] && !gMC->IsNewTrack()) {
323       TList *histnew = dynamic_cast<TList*>((*fQAList)[id]);
324       hist = dynamic_cast<TH1F*>(histnew->FindObject("hEnIn"));
325       hist->Fill(TMath::Log10(energy));
326       hist = dynamic_cast<TH1F*>(histnew->FindObject("hZIn"));
327       hist->Fill(x[2]);
328     }
329     fOldId=id;
330   }
331 }
332
333 //_______________________________________________________________________
334 void AliMCQA::AddModuleName()
335 {
336   //
337   // Add function DrawModuleName to the module frequency histogram pad
338   //
339   static TVirtualPad* oldpad=0;
340   if(oldpad!=gPad) {
341     gPad->GetCanvas()->FeedbackMode(kTRUE);
342     if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawModuleName()");
343     DrawPaveLabel(fMPaveLabel);
344     oldpad=gPad;
345   }
346 }
347
348 //_______________________________________________________________________
349 void AliMCQA::AddVolumeName()
350 {
351   //
352   // Add function DrawVolumeName to the volume frequency histogram pad
353   //
354   static TVirtualPad* oldpad=0;
355   if(oldpad!=gPad) {
356     gPad->GetCanvas()->FeedbackMode(kTRUE);
357     if(gPad) gPad->AddExec("ex","gAlice->GetMCQA()->DrawVolumeName()");
358     DrawPaveLabel(fVPaveLabel);
359     oldpad=gPad;
360   }
361 }
362
363 //_______________________________________________________________________
364 void AliMCQA::DrawPaveLabel(TPaveLabel *&pv)
365 {
366   //
367   // Draws the PaveLabel with the meaning of the bin
368   //
369   float uxmin = gPad->GetUxmin();
370   float uxmax = gPad->GetUxmax();
371   float uymin = gPad->GetUymin();
372   float uymax = gPad->GetUymax();
373   float lx = uxmax-uxmin;
374   float ly = uymax-uymin;
375
376   if(pv) delete pv;
377   pv 
378     = new TPaveLabel(uxmin+0.05*lx,uymax-(0.05+0.1)*ly,
379                      uxmin+(0.05+0.2)*lx,uymax-0.05*ly,
380                      "");
381   pv->Draw();
382 }
383
384 //_______________________________________________________________________
385 Int_t AliMCQA::GetHBin(const char* hname)
386 {
387   //
388   // Get the bin where the cursor is
389   //
390   TList *dir = gDirectory->GetList();
391   TH1 *h=dynamic_cast<TH1*>(dir->FindObject(hname));
392   
393
394   int px = gPad->GetEventX();
395   Float_t upx = gPad->AbsPixeltoX(px);
396   Float_t x = gPad->PadtoX(upx);
397     
398   return h->GetXaxis()->FindBin(x);
399 }
400
401 //_______________________________________________________________________
402 void AliMCQA::DrawModuleName()
403 {
404   //
405   // Writes the name of the module of the bin where the cursor is
406   //
407   TObject *select = gPad->GetSelected();
408   if(!select) return;
409   
410   Int_t binx = GetHBin("hMCMcalls");
411   if(0<binx && binx<=fNdets) {
412     char lab[15];
413     strcpy(lab,dynamic_cast<TNamed*>((*fModNames)[binx-1])->GetName());
414     fMPaveLabel->SetLabel(lab);
415   
416     gPad->Modified();
417     gPad->Update();
418   }
419 }
420
421 //_______________________________________________________________________
422 void AliMCQA::DrawVolumeName()
423 {
424   //
425   // Writes the name of the volume:module of the bin where the cursor is
426   //
427   TObject *select = gPad->GetSelected();
428   if(!select) return;
429
430   Int_t binx = GetHBin("hMCVcalls");
431   if(0<binx && binx<=fNvolumes) {
432     char lab[20];
433     sprintf(lab,"%s: %s",dynamic_cast<TNamed*>((*fVolNames)[binx-1])->GetName(),
434             dynamic_cast<TNamed*>((*fVolNames)[binx-1])->GetTitle());
435     fVPaveLabel->SetLabel(lab);
436     
437     gPad->Modified();
438     gPad->Update();
439   }
440 }
441
442