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