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