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