AliRoot Core  ee782a0 (ee782a0)
AliTransportMonitor.cxx
Go to the documentation of this file.
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 // Class that can be plugged in the simulation to monitor transport timing per
17 // particle for each geometry volume.
18 //
19 // andrei.gheata@cern.ch
20 
21 #include "AliTransportMonitor.h"
22 
23 #include "TDatabasePDG.h"
24 #include "TCanvas.h"
25 #include "TMath.h"
26 #include "TFile.h"
27 #include "TROOT.h"
28 #include "THashList.h"
29 #include "TH2F.h"
30 #include "TGeoManager.h"
31 #include "AliPDG.h"
32 #include "TVirtualMC.h"
33 #include <RVersion.h>
34 
35 ClassImp(AliTransportMonitor)
36 ClassImp(AliTransportMonitor::AliTransportMonitorVol)
37 ClassImp(AliTransportMonitor::AliTransportMonitorVol::AliPMonData)
38 
39 typedef AliTransportMonitor::AliTransportMonitorVol::AliPMonData PMonData;
40 
41 //______________________________________________________________________________
43  :TNamed(),
44  fNtypes(0),
45  fTotalTime(0),
46  fNSteps(0),
47  fPData(0),
48  fTimeRZ(0),
49  fParticles()
50 {
51 // Default constructor
52 }
53 
54 //______________________________________________________________________________
56 {
57 // Destructor
58  delete [] fPData;
59  delete fTimeRZ;
60 }
61 
62 //______________________________________________________________________________
64  Int_t pdg,
65  Double_t energy,
66  Double_t dt,
67  Double_t x, Double_t y, Double_t z)
68 {
69 // This method is called at each N steps to store timing info.
70  PMonData &data = GetPMonData(pdg);
71  data.fEdt += energy*dt;
72  data.fTime += dt;
73  fTotalTime += dt;
74  fNSteps += 1.;
75  if (!fTimeRZ) {
76  Bool_t status = TH1::AddDirectoryStatus();
77  TH1::AddDirectory(kFALSE);
78  fTimeRZ = new TH2F("h2rz", "", 100, -5000., 5000, 100, 0., 1000.);
79  TH1::AddDirectory(status);
80  }
81  Double_t r = TMath::Sqrt(x*x+y*y);
82  fTimeRZ->Fill(z,r,dt);
83 }
84 
85 //______________________________________________________________________________
87 {
88 // Retrieve stored monitoring object for a given pdg type. If not existing
89 // create one.
90 
91  // The object could have been retrieved from file, in which case we have to
92  // build the map.
93  //
94  // unknown heavy fragment ?
95  // TParticlePDG* pdgP = (TDatabasePDG::Instance())->GetParticle(pdg);
96  Int_t apdg = TMath::Abs(pdg);
97  if ((apdg > 10000)
98  && (apdg != 1000010020)
99  && (apdg != 1000010030)
100  && (apdg != 1000020030)
101  && (apdg != 1000020040)
102  && (apdg != 50000050)
103  && (apdg != 50000051)
104  )
105  pdg = 1111111111;
106 
107  PMonData *data;
108  if (fNtypes) {
109  if (fParticles.empty()) {
110  for (Int_t i=0; i<fNtypes; i++) {
111  data = &fPData[i];
112  fParticles[i] = data->fPDG;
113  }
114  }
115  ParticleMapIt_t it = fParticles.find(pdg);
116  if (it == fParticles.end()) {
117  data = &fPData[fNtypes];
118  data->fPDG = pdg;
119  fParticles[pdg] = fNtypes++;
120  } else {
121  data = &fPData[it->second];
122  }
123  } else {
124  TDatabasePDG *pdgDB = TDatabasePDG::Instance();
125  if (!pdgDB->ParticleList()) AliPDG::AddParticlesToPdgDataBase();
126  Int_t size = pdgDB->ParticleList()->GetSize();
127  // account for heavy fragments coded as "1111111111"
128  //
129  fPData = new PMonData[size+10];
130  data = &fPData[fNtypes];
131  data->fPDG = pdg;
132  fParticles[pdg] = fNtypes++;
133  }
134  return *data;
135 }
136 
138 {
139  //
140  // Merging
141  //
142  fTotalTime = (fTotalTime + volM->GetTotalTime());
143  if (fTimeRZ && volM->GetHistogram()) {
144  fTimeRZ->Add(volM->GetHistogram());
145  } else if (volM->GetHistogram()) {
146  fTimeRZ = (TH2F*)(volM->GetHistogram()->Clone());
147  }
148 
149  Int_t ntypes = volM->GetNtypes();
150  for (Int_t i = 0; i < ntypes; i++) {
151  Int_t pdg = volM->GetPDG(i);
152  PMonData &data = GetPMonData(pdg);
153  data.fEdt += (volM->GetEmed(i) * volM->GetTotalTime());
154  data.fTime += (volM->GetTime(i));
155  }
156 }
157 
158 ClassImp(AliTransportMonitor)
159 
160 //______________________________________________________________________________
162  :TObject(),
163  fTotalTime(0),
164  fTimer(),
165  fVolumeMon(0)
166 {
167 // Default constructor
168 }
169 
170 //______________________________________________________________________________
172  :TObject(),
173  fTotalTime(0),
174  fTimer(),
175  fVolumeMon(0)
176 {
177 // Default constructor
178  fVolumeMon = new TObjArray(nvolumes);
179  fVolumeMon->SetOwner();
180  for (Int_t i=0; i<nvolumes; i++) {
182  if (TVirtualMC::GetMC()) volMon->SetName(TVirtualMC::GetMC()->VolName(i));
183  fVolumeMon->Add(volMon);
184  }
185 }
186 
187 //______________________________________________________________________________
189 {
190 // Destructor
191  delete fVolumeMon;
192 }
193 
194 //______________________________________________________________________________
195 void AliTransportMonitor::Print(Option_t *volName) const
196 {
197 // Inspect the timing statistics for a single volume or for all the setup
198  Int_t uid = -1;
199  Int_t ntotal = 0;
200  if (!fVolumeMon || !(ntotal=fVolumeMon->GetEntriesFast())) {
201  Info("Inspect", "Transport monitor is empty !");
202  return;
203  }
204  if (strlen(volName)) {
205  TString svname = volName;
206  Int_t i = 0;
207  for (i=1; i<ntotal; i++) if (svname == fVolumeMon->At(i)->GetName()) break;
208  if (i==ntotal) {
209  Error("Inspect", "No monitoring info stored for volume %s", volName);
210  return;
211  }
212  uid = i;
214  Int_t ntypes = volMon->GetNtypes();
215  if (!ntypes) {
216  Info("Inspect", "No particles crossed volume %s", volName);
217  return;
218  }
219  Double_t *timeperpart = new Double_t[ntypes];
220  Int_t *isort = new Int_t[ntypes];
221  Double_t timepervol = 0.;
222  for (i=0; i<ntypes; i++) {
223  timeperpart[i] = volMon->GetTime(i);
224  timepervol += timeperpart[i];
225  }
226  printf("Volume %s: Transport time: %g%% of %g %g [s]\n", volMon->GetName(), 100.*timepervol/fTotalTime, fTotalTime,
227  volMon->GetTotalTime());
228  TMath::Sort(ntypes, timeperpart, isort, kTRUE);
229  TString particle;
230  TDatabasePDG *pdgDB = TDatabasePDG::Instance();
231  if (!pdgDB->ParticleList()) AliPDG::AddParticlesToPdgDataBase();
232  for (i=0; i<ntypes; i++) {
233  timeperpart[i] /= timepervol;
234  TParticlePDG* pdgP = pdgDB->GetParticle(volMon->GetPDG(isort[i]));
235  if (pdgP) {
236  particle = pdgDB->GetParticle(volMon->GetPDG(isort[i]))->GetName();
237  } else {
238  particle = Form("pdg code not in DB: %d", volMon->GetPDG(isort[i]));
239  }
240  printf(" %s: %g%% mean energy: %g\n", particle.Data(), 100.*timeperpart[i], volMon->GetEmed(i));
241  }
242  if (volMon->GetHistogram()) {
243  TCanvas *c1 = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("crz");
244  if (!c1) c1 = new TCanvas("crz");
245  c1->cd();
246  volMon->GetHistogram()->GetXaxis()->SetTitle("z [cm]");
247  volMon->GetHistogram()->GetYaxis()->SetTitle("r [cm]");
248  volMon->GetHistogram()->SetTitle(Form("RZ plot weighted by time spent in %s",volMon->GetName()));
249  volMon->GetHistogram()->Draw();
250  }
251  return;
252  }
253  // General view
254  TIter next(fVolumeMon);
255  AliTransportMonitorVol *volMon;
256  Int_t ncrossed = 0;
257  TH1F *hnames = new TH1F("volume_timing", "relative volume timing", 3,0,3);
258  hnames->SetStats(0);
259  hnames->SetFillColor(38);
260  #if ROOT_VERSION_CODE < ROOT_VERSION(6,4,0)
261  hnames->SetBit(TH1::kCanRebin);
262  #endif
263  while ((volMon=(AliTransportMonitorVol*)next())) {
264  if (volMon->GetNtypes()) {
265  hnames->Fill(volMon->GetName(), volMon->GetTotalTime());
266  ncrossed++;
267  }
268  }
269 
270  hnames->LabelsDeflate();
271  hnames->GetXaxis()->LabelsOption(">");
272 
273  TCanvas *c = (TCanvas*)gROOT->GetListOfCanvases()->FindObject("cvol_timing");
274  if (!c) c = new TCanvas("cvol_timing");
275  c->cd();
276  c->SetLogy();
277  hnames->Draw();
278 
279  printf("=============================================================================\n");
280  printf("Effective transport time: %6.2f minutes\n", fTotalTime/60.);
281  printf("Number of crossed volumes: %d from %d\n", ncrossed, fVolumeMon->GetEntriesFast());
282  printf("=============================================================================\n");
283 }
284 
285 //______________________________________________________________________________
287 {
288 // Reset timer for zero-length steps
289  fTimer.Stop();
290  fTimer.Start(kTRUE);
291 }
292 
293 //______________________________________________________________________________
295  Int_t pdg,
296  Double_t energy,
297  Double_t x, Double_t y, Double_t z)
298 {
299 // This method is called at each N steps to store timing info.
300  fTimer.Stop();
301  Double_t dt = fTimer.RealTime();
302  fTotalTime += dt;
304  volMon->StepInfo(pdg,energy,dt,x,y,z);
305  fTimer.Start(kTRUE);
306 }
307 
308 //______________________________________________________________________________
310 {
311 // Start collecting timing information
312  if (fTotalTime > 0) {
313  Info("Start", "Cannot start twice");
314  return;
315  }
316  fTimer.Start(kTRUE);
317 }
318 
319 //______________________________________________________________________________
321 {
322 // Stop the global timer
323  fTimer.Stop();
324 }
325 
326 //______________________________________________________________________________
328 {
329 // Export information to file
330  TFile *file = TFile::Open(fname, "RECREATE");
331  Write();
332  file->Write();
333  file->Close();
334 }
335 //______________________________________________________________________________
337 {
338 
339  // merge with monitor
340  if (!fVolumeMon)
341  {
342  TObjArray* arr = mergeMon->GetVolumes();
343  Int_t nvol = arr->GetEntriesFast();
344  fVolumeMon = new TObjArray(nvol);
345  fVolumeMon->SetOwner();
346  for (Int_t i = 0; i < nvol; i++) {
348  volMon->SetName(arr->At(i)->GetName());
349  fVolumeMon->Add(volMon);
350  }
351  } // first time
352 
353 
354  Int_t n = fVolumeMon->GetEntriesFast();
355  TObjArray* mergeVols = mergeMon->GetVolumes();
356  fTotalTime = 0;
357  for (Int_t i = 0; i < n; i++)
358  {
360  AliTransportMonitorVol *volMon2 = (AliTransportMonitorVol*)mergeVols->At(i);
361  volMon1->Merge(volMon2);
362  fTotalTime += (volMon1->GetTotalTime());
363  }
364 }
365 //______________________________________________________________________________
367 {
368 // Import information from a file
369  TFile *file = TFile::Open(fname);
370  if (!file) {
371  ::Error("Import", "File %s could not be opened", fname);
372  return 0;
373  }
374  AliTransportMonitor *mon = (AliTransportMonitor *)file->Get("AliTransportMonitor");
375  if (!mon) {
376  ::Error("Import", "No AliTransportMonitor object found n file %s", fname);
377  return 0;
378  }
380  return mon;
381 }
382 
383 
TDatime dt
Definition: pdc06_config.C:130
printf("Chi2/npoints = %f\n", TMath::Sqrt(chi2/npoints))
TFile * Open(const char *filename, Long64_t &nevents)
TObjArray * GetVolumes() const
#define TObjArray
void StepInfo(Int_t pdg, Double_t energy, Double_t dt, Double_t x, Double_t y, Double_t z)
void Merge(AliTransportMonitor *mergeMon)
TROOT * gROOT
static void AddParticlesToPdgDataBase()
Definition: AliPDG.cxx:31
void StepInfo(Int_t volId, Int_t pdg, Double_t energy, Double_t x, Double_t y, Double_t z)
void Print(Option_t *volName="") const
TObjArray * fVolumeMon
Global timer.
AliPMonData * GetPMonData(Int_t itype) const
static AliTransportMonitor * Import(const char *fname)
void Export(const char *fname)
char * fname
void Merge(AliTransportMonitorVol *volM)