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