added vertex reco systematic stuff.
[u/mrichter/AliRoot.git] / PWG0 / dNdEta / AlidNdEtaSystematicsSelector.cxx
1 /* $Id$ */
2
3 #include "AlidNdEtaSystematicsSelector.h"
4
5 #include <TStyle.h>
6 #include <TSystem.h>
7 #include <TCanvas.h>
8 #include <TVector3.h>
9 #include <TChain.h>
10 #include <TFile.h>
11 #include <TH2F.h>
12 #include <TH1F.h>
13 #include <TParticle.h>
14
15 #include <AliLog.h>
16 #include <AliESD.h>
17 #include <AliStack.h>
18 #include <AliHeader.h>
19 #include <AliGenEventHeader.h>
20 #include <AliGenPythiaEventHeader.h>
21 #include <AliGenCocktailEventHeader.h>
22
23
24 #include "esdTrackCuts/AliESDtrackCuts.h"
25 #include "AliPWG0Helper.h"
26 #include "dNdEta/AlidNdEtaCorrection.h"
27
28 ClassImp(AlidNdEtaSystematicsSelector)
29
30 AlidNdEtaSystematicsSelector::AlidNdEtaSystematicsSelector() :
31   AliSelectorRL(),
32   fSecondaries(0),
33   fSigmaVertex(0),
34   fEsdTrackCuts(0),
35   fPIDParticles(0),
36   fPIDTracks(0)
37 {
38   //
39   // Constructor. Initialization of pointers
40   //
41
42   for (Int_t i=0; i<4; ++i)
43     fdNdEtaCorrectionSpecies[i] = 0;
44
45   for (Int_t i=0; i<3; ++i)
46     fdNdEtaCorrectionVertexReco[i] = 0;
47
48 }
49
50 AlidNdEtaSystematicsSelector::~AlidNdEtaSystematicsSelector()
51 {
52   //
53   // Destructor
54   //
55 }
56
57 void AlidNdEtaSystematicsSelector::Begin(TTree* tree)
58 {
59   // Begin function
60
61   AliSelectorRL::Begin(tree);
62
63   ReadUserObjects(tree);
64 }
65
66 void AlidNdEtaSystematicsSelector::ReadUserObjects(TTree* tree)
67 {
68   // read the user objects, called from slavebegin and begin
69
70   if (!fEsdTrackCuts && fInput)
71     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (fInput->FindObject("AliESDtrackCuts"));
72
73   if (!fEsdTrackCuts && tree)
74     fEsdTrackCuts = dynamic_cast<AliESDtrackCuts*> (tree->GetUserInfo()->FindObject("AliESDtrackCuts"));
75
76   if (!fEsdTrackCuts)
77      AliDebug(AliLog::kError, "ERROR: Could not read EsdTrackCuts from input list.");
78 }
79
80 void AlidNdEtaSystematicsSelector::SlaveBegin(TTree* tree)
81 {
82   // The SlaveBegin() function is called after the Begin() function.
83   // When running with PROOF SlaveBegin() is called on each slave server.
84   // The tree argument is deprecated (on PROOF 0 is passed).
85
86   AliSelector::SlaveBegin(tree);
87
88   ReadUserObjects(tree);
89
90   TString option(GetOption());
91
92   printf("Running AlidNdEtaSystematicsSelector with options %s\n", option.Data());
93
94   if (option.Contains("secondaries"))
95   {
96     fSecondaries = new TH2F("fSecondaries", "fSecondaries;Case;N", 8, 0.5, 8.5, 2001, -0.5, 2000.5);
97     fSecondaries->GetXaxis()->SetBinLabel(1, "All Primaries");
98     fSecondaries->GetXaxis()->SetBinLabel(2, "All Secondaries");
99     fSecondaries->GetXaxis()->SetBinLabel(3, "Primaries pT > 0.3 GeV/c");
100     fSecondaries->GetXaxis()->SetBinLabel(4, "Secondaries pT > 0.3 GeV/c");
101     fSecondaries->GetXaxis()->SetBinLabel(5, "Tracks from Primaries");
102     fSecondaries->GetXaxis()->SetBinLabel(6, "Tracks from Secondaries");
103     fSecondaries->GetXaxis()->SetBinLabel(7, "Accepted Tracks from Primaries");
104     fSecondaries->GetXaxis()->SetBinLabel(8, "Accepted Tracks from Secondaries");
105   }
106
107   if (option.Contains("particle-composition"))
108   {
109     for (Int_t i=0; i<4; ++i)
110     {
111       TString name;
112       name.Form("correction_%d", i);
113       fdNdEtaCorrectionSpecies[i] = new AlidNdEtaCorrection(name, name);
114     }
115   }
116
117   if (option.Contains("sigma-vertex"))
118   {
119     fSigmaVertex = new TH1F("fSigmaVertex", "fSigmaVertex;Nsigma2vertex;NacceptedTracks", 50, 0.05, 5.05);
120     printf("WARNING: sigma-vertex analysis enabled. This will produce weird results in the AliESDtrackCuts histograms\n");
121   }
122
123   if (option.Contains("vertexreco")) {
124     fdNdEtaCorrectionVertexReco[0] = new AlidNdEtaCorrection("vertexRecoND", "vertexRecoND");
125     fdNdEtaCorrectionVertexReco[1] = new AlidNdEtaCorrection("vertexRecoSD", "vertexRecoSD");
126     fdNdEtaCorrectionVertexReco[2] = new AlidNdEtaCorrection("vertexRecoDD", "vertexRecoDD");
127   }
128
129   fPIDParticles = new TH1F("pid_particles", "PID of generated primary particles", 10001, -5000.5, 5000.5);
130   fPIDTracks = new TH1F("pid_tracks", "MC PID of reconstructed tracks", 10001, -5000.5, 5000.5);
131 }
132
133 Bool_t AlidNdEtaSystematicsSelector::Process(Long64_t entry)
134 {
135   // The Process() function is called for each entry in the tree (or possibly
136   // keyed object in the case of PROOF) to be processed. The entry argument
137   // specifies which entry in the currently loaded tree is to be processed.
138   // It can be passed to either TTree::GetEntry() or TBranch::GetEntry()
139   // to read either all or the required parts of the data. When processing
140   // keyed objects with PROOF, the object is already loaded and is available
141   // via the fObject pointer.
142   //
143   // This function should contain the "body" of the analysis. It can contain
144   // simple or elaborate selection criteria, run algorithms on the data
145   // of the event and typically fill histograms.
146
147   // WARNING when a selector is used with a TChain, you must use
148   //  the pointer to the current TTree to call GetEntry(entry).
149   //  The entry is always the local entry number in the current tree.
150   //  Assuming that fTree is the pointer to the TChain being processed,
151   //  use fTree->GetTree()->GetEntry(entry).
152
153   if (AliSelectorRL::Process(entry) == kFALSE)
154     return kFALSE;
155
156   // Check prerequisites
157   if (!fESD)
158   {
159     AliDebug(AliLog::kError, "ESD branch not available");
160     return kFALSE;
161   }
162
163   if (!fEsdTrackCuts)
164   {
165     AliDebug(AliLog::kError, "fESDTrackCuts not available");
166     return kFALSE;
167   }
168
169   AliStack* stack = GetStack();
170   if (!stack)
171   {
172     AliDebug(AliLog::kError, "Stack not available");
173     return kFALSE;
174   }
175
176   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
177
178   if (fdNdEtaCorrectionSpecies[0])
179     FillCorrectionMaps(list);
180
181   if (fSecondaries)
182     FillSecondaries();
183
184   if (fSigmaVertex)
185     FillSigmaVertex();
186
187   // stuff for vertex reconstruction correction systematics
188   if (fdNdEtaCorrectionVertexReco[0]) {
189     AliHeader* header = GetHeader();
190     if (!header) {
191       AliDebug(AliLog::kError, "Header not available");
192       return kFALSE;
193     }
194
195     // is there a smart way to check if it is really a pythia header?
196     AliGenCocktailEventHeader* genHeaders = dynamic_cast<AliGenCocktailEventHeader*>(header->GenEventHeader());
197     if (!genHeaders)
198       return kTRUE;
199
200     TList* headerList = genHeaders->GetHeaders();
201     if (!headerList)
202       return kTRUE;
203     
204     AliGenPythiaEventHeader* genHeader = 0;
205     for (Int_t i=0; i<headerList->GetEntries(); i++) {
206       genHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
207       if (genHeader)
208         break;
209     }        
210     if (!genHeader)
211       return kTRUE;
212     
213     Int_t processtype = genHeader->ProcessType();      
214
215     // get the MC vertex
216     TArrayF vtxMC(3);
217     genHeader->PrimaryVertex(vtxMC);
218
219     Int_t nGoodTracks = list->GetEntries();
220
221     Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD);
222     Bool_t vertexReconstructed = AliPWG0Helper::IsVertexReconstructed(fESD);
223
224     // genHeader->Print();
225     // printf("  heps %d %d %d %f \n", eventTriggered, vertexReconstructed, processtype, vtxMC[2]);
226
227     // non diffractive
228     if (processtype!=92 && processtype!=93 && processtype!=94) { 
229       fdNdEtaCorrectionVertexReco[0]->FillEvent(vtxMC[2], nGoodTracks);
230       if (eventTriggered) {
231         fdNdEtaCorrectionVertexReco[0]->FillEventWithTrigger(vtxMC[2], nGoodTracks);
232         if (vertexReconstructed)
233           fdNdEtaCorrectionVertexReco[0]->FillEventWithTriggerWithReconstructedVertex(vtxMC[2], nGoodTracks);
234       }
235     }
236
237     // single diffractive
238     if (processtype==92 || processtype==93) { 
239       fdNdEtaCorrectionVertexReco[1]->FillEvent(vtxMC[2], nGoodTracks);
240       if (eventTriggered) {
241         fdNdEtaCorrectionVertexReco[1]->FillEventWithTrigger(vtxMC[2], nGoodTracks);
242         if (vertexReconstructed)
243           fdNdEtaCorrectionVertexReco[1]->FillEventWithTriggerWithReconstructedVertex(vtxMC[2], nGoodTracks);
244       }
245     }
246
247     // double diffractive
248     if (processtype==94) { 
249       fdNdEtaCorrectionVertexReco[2]->FillEvent(vtxMC[2], nGoodTracks);
250       if (eventTriggered) {
251         fdNdEtaCorrectionVertexReco[2]->FillEventWithTrigger(vtxMC[2], nGoodTracks);
252         if (vertexReconstructed)
253           fdNdEtaCorrectionVertexReco[2]->FillEventWithTriggerWithReconstructedVertex(vtxMC[2], nGoodTracks);
254       }
255     }
256   }
257
258   delete list;
259   list = 0;
260
261   return kTRUE;
262 }
263
264 void AlidNdEtaSystematicsSelector::FillCorrectionMaps(TObjArray* listOfTracks)
265 {
266   // fills the correction maps for different particle species
267
268   AliStack* stack = GetStack();
269   AliHeader* header = GetHeader();
270
271   Bool_t vertexReconstructed = AliPWG0Helper::IsVertexReconstructed(fESD);
272   Bool_t eventTriggered = AliPWG0Helper::IsEventTriggered(fESD);
273
274   // get the MC vertex
275   AliGenEventHeader* genHeader = header->GenEventHeader();
276
277   TArrayF vtxMC(3);
278   genHeader->PrimaryVertex(vtxMC);
279
280   // loop over mc particles
281   Int_t nPrim  = stack->GetNprimary();
282
283   for (Int_t iMc = 0; iMc < nPrim; ++iMc)
284   {
285     TParticle* particle = stack->Particle(iMc);
286
287     if (!particle)
288     {
289       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (mc loop).", iMc));
290       continue;
291     }
292
293     if (AliPWG0Helper::IsPrimaryCharged(particle, nPrim) == kFALSE)
294       continue;
295
296     Float_t eta = particle->Eta();
297     Float_t pt = particle->Pt();
298
299     Int_t id = -1;
300     switch (TMath::Abs(particle->GetPdgCode()))
301     {
302       case 211: id = 0; break;
303       case 321: id = 1; break;
304       case 2212: id = 2; break;
305       case 11:
306         {
307           /*if (pt < 0.1 && particle->GetMother(0) > -1)
308           {
309             TParticle* mother = stack->Particle(particle->GetMother(0));
310             printf("Mother pdg code is %d\n", mother->GetPdgCode());
311           } */
312           //particle->Dump();
313           //if (particle->GetUniqueID() == 1)
314
315           //printf("Found illegal particle mc id = %d file = %s event = %d\n", iMc,  fTree->GetCurrentFile()->GetName(), fTree->GetTree()->GetReadEntry());
316         }
317       default: id = 3; break;
318     }
319
320     if (vertexReconstructed)
321     {
322       fdNdEtaCorrectionSpecies[id]->FillParticle(vtxMC[2], eta, pt);
323       //if (pt < 0.1)
324         fPIDParticles->Fill(particle->GetPdgCode());
325     }
326
327     fdNdEtaCorrectionSpecies[id]->FillParticleAllEvents(eta, pt);
328     if (eventTriggered)
329       fdNdEtaCorrectionSpecies[id]->FillParticleWhenEventTriggered(eta, pt);
330   }// end of mc particle
331
332   // loop over esd tracks
333   TIterator* iter = listOfTracks->MakeIterator();
334   TObject* obj = 0;
335   while ((obj = iter->Next()) != 0)
336   {
337     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
338     if (!esdTrack)
339       continue;
340
341     // using the properties of the mc particle
342     Int_t label = TMath::Abs(esdTrack->GetLabel());
343     TParticle* particle = stack->Particle(label);
344     if (!particle)
345     {
346       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
347       continue;
348     }
349
350     TParticle* mother = particle;
351     // find primary particle that created this particle
352     while (label >= nPrim)
353     {
354       //printf("Particle %d (pdg %d) is not a primary. Let's check its mother %d\n", label, mother->GetPdgCode(), mother->GetMother(0));
355
356       if (mother->GetMother(0) == -1)
357       {
358         AliDebug(AliLog::kError, Form("UNEXPECTED: Could not find mother of secondary particle %d.", label));
359         mother = 0;
360         break;
361       }
362
363       label = mother->GetMother(0);
364
365       mother = stack->Particle(label);
366       if (!mother)
367       {
368         AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (find mother loop).", label));
369         break;
370       }
371     }
372
373     if (!mother)
374       continue;
375
376     //printf("We continue with particle %d (pdg %d)\n", label, mother->GetPdgCode());
377
378     Int_t id = -1;
379     switch (TMath::Abs(mother->GetPdgCode()))
380     {
381       case 211:  id = 0; break;
382       case 321:  id = 1; break;
383       case 2212: id = 2; break;
384       default:   id = 3; break;
385     }
386
387     if (vertexReconstructed)
388     {
389       fdNdEtaCorrectionSpecies[id]->FillParticleWhenMeasuredTrack(vtxMC[2], particle->Eta(), particle->Pt());
390       //if (particle->Pt() < 0.1)
391         fPIDTracks->Fill(particle->GetPdgCode());
392     }
393   } // end of track loop
394
395   Int_t nGoodTracks = listOfTracks->GetEntries();
396
397   for (Int_t i=0; i<4; ++i)
398   {
399     fdNdEtaCorrectionSpecies[i]->FillEvent(vtxMC[2], nGoodTracks);
400     if (eventTriggered)
401     {
402       fdNdEtaCorrectionSpecies[i]->FillEventWithTrigger(vtxMC[2], nGoodTracks);
403       if (vertexReconstructed)
404         fdNdEtaCorrectionSpecies[i]->FillEventWithTriggerWithReconstructedVertex(vtxMC[2], nGoodTracks);
405     }
406   }
407
408   delete iter;
409   iter = 0;
410 }
411
412 void AlidNdEtaSystematicsSelector::FillSecondaries()
413 {
414   // fills the secondary histograms
415
416   AliStack* stack = GetStack();
417
418   Int_t particlesPrimaries = 0;
419   Int_t particlesSecondaries = 0;
420   Int_t particlesPrimariesPtCut = 0;
421   Int_t particlesSecondariesPtCut = 0;
422
423   for (Int_t iParticle = 0; iParticle < stack->GetNtrack(); iParticle++)
424   {
425     TParticle* particle = stack->Particle(iParticle);
426     if (!particle)
427     {
428       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (particle loop).", iParticle));
429       continue;
430     }
431
432     if (TMath::Abs(particle->Eta()) > 0.9)
433       continue;
434
435     Bool_t isPrimary = kFALSE;
436     if (iParticle < stack->GetNprimary())
437     {
438       if (AliPWG0Helper::IsPrimaryCharged(particle, stack->GetNprimary()) == kFALSE)
439         continue;
440
441       isPrimary = kTRUE;
442     }
443
444     if (isPrimary)
445       particlesPrimaries++;
446     else
447       particlesSecondaries++;
448
449     if (particle->Pt() < 0.3)
450       continue;
451
452     if (isPrimary)
453       particlesPrimariesPtCut++;
454     else
455       particlesSecondariesPtCut++;
456   }
457
458   fSecondaries->Fill(1, particlesPrimaries);
459   fSecondaries->Fill(2, particlesSecondaries);
460   fSecondaries->Fill(3, particlesPrimariesPtCut);
461   fSecondaries->Fill(4, particlesSecondariesPtCut);
462
463   Int_t allTracksPrimaries = 0;
464   Int_t allTracksSecondaries = 0;
465   Int_t acceptedTracksPrimaries = 0;
466   Int_t acceptedTracksSecondaries = 0;
467
468   for (Int_t iTrack = 0; iTrack < fESD->GetNumberOfTracks(); iTrack++)
469   {
470     AliESDtrack* esdTrack = fESD->GetTrack(iTrack);
471     if (!esdTrack)
472     {
473       AliDebug(AliLog::kError, Form("UNEXPECTED: Could not get track %d.", iTrack));
474       continue;
475     }
476
477     Int_t label = TMath::Abs(esdTrack->GetLabel());
478     TParticle* particle = stack->Particle(label);
479     if (!particle)
480     {
481       AliDebug(AliLog::kError, Form("UNEXPECTED: particle with label %d not found in stack (track loop).", label));
482       continue;
483     }
484
485    if (label < stack->GetNprimary())
486       allTracksPrimaries++;
487     else
488       allTracksSecondaries++;
489
490     if (!fEsdTrackCuts->AcceptTrack(esdTrack))
491       continue;
492
493     if (label < stack->GetNprimary())
494       acceptedTracksPrimaries++;
495     else
496       acceptedTracksSecondaries++;
497   }
498
499   fSecondaries->Fill(5, allTracksPrimaries);
500   fSecondaries->Fill(6, allTracksSecondaries);
501   fSecondaries->Fill(7, acceptedTracksPrimaries);
502   fSecondaries->Fill(8, acceptedTracksSecondaries);
503
504   //printf("P = %d, S = %d, P_pt = %d, S_pt = %d, T_P = %d, T_S = %d, T_P_acc = %d, T_S_acc = %d\n", particlesPrimaries, particlesSecondaries, particlesPrimariesPtCut, particlesSecondariesPtCut, allTracksPrimaries, allTracksSecondaries, acceptedTracksPrimaries, acceptedTracksSecondaries);
505 }
506
507 void AlidNdEtaSystematicsSelector::FillSigmaVertex()
508 {
509   // fills the fSigmaVertex histogram
510
511   // save the old value
512   Float_t oldSigmaVertex = fEsdTrackCuts->GetMinNsigmaToVertex();
513
514   // set to maximum
515   fEsdTrackCuts->SetMinNsigmaToVertex(5);
516
517   TObjArray* list = fEsdTrackCuts->GetAcceptedTracks(fESD);
518
519   TIterator* iter = list->MakeIterator();
520   TObject* obj = 0;
521   while ((obj = iter->Next()) != 0)
522   {
523     AliESDtrack* esdTrack = dynamic_cast<AliESDtrack*> (obj);
524     if (!esdTrack)
525       continue;
526
527     Float_t sigma2Vertex = fEsdTrackCuts->GetSigmaToVertex(esdTrack);
528
529     for (Double_t nSigma = 0.1; nSigma < 5.05; nSigma += 0.1)
530     {
531       if (sigma2Vertex < nSigma)
532         fSigmaVertex->Fill(nSigma);
533     }
534   }
535
536   delete iter;
537   iter = 0;
538
539   delete list;
540   list = 0;
541
542   // set back the old value
543   fEsdTrackCuts->SetMinNsigmaToVertex(oldSigmaVertex);
544 }
545
546 void AlidNdEtaSystematicsSelector::SlaveTerminate()
547 {
548   // The SlaveTerminate() function is called after all entries or objects
549   // have been processed. When running with PROOF SlaveTerminate() is called
550   // on each slave server.
551
552   AliSelector::SlaveTerminate();
553
554   // Add the histograms to the output on each slave server
555   if (!fOutput)
556   {
557     AliDebug(AliLog::kError, Form("ERROR: Output list not initialized."));
558     return;
559   }
560
561   if (fSecondaries)
562     fOutput->Add(fSecondaries);
563
564   for (Int_t i=0; i<4; ++i)
565     if (fdNdEtaCorrectionSpecies[i])
566       fOutput->Add(fdNdEtaCorrectionSpecies[i]);
567
568   if (fSigmaVertex)
569     fOutput->Add(fSigmaVertex);
570
571   for (Int_t i=0; i<3; ++i)
572     if (fdNdEtaCorrectionVertexReco[i])
573       fOutput->Add(fdNdEtaCorrectionVertexReco[i]);
574
575 }
576
577 void AlidNdEtaSystematicsSelector::Terminate()
578 {
579   // The Terminate() function is the last function to be called during
580   // a query. It always runs on the client, it can be used to present
581   // the results graphically or save the results to file.
582
583   AliSelector::Terminate();
584
585   fSecondaries = dynamic_cast<TH2F*> (fOutput->FindObject("fSecondaries"));
586   for (Int_t i=0; i<4; ++i)
587     fdNdEtaCorrectionSpecies[i] = dynamic_cast<AlidNdEtaCorrection*> (fOutput->FindObject(Form("correction_%d", i)));
588   fSigmaVertex = dynamic_cast<TH1F*> (fOutput->FindObject("fSigmaVertex"));
589
590   if (fPIDParticles)
591   {
592     TDatabasePDG* pdgDB = new TDatabasePDG;
593
594     for (Int_t i=0; i <= fPIDParticles->GetNbinsX()+1; ++i)
595       if (fPIDParticles->GetBinContent(i) > 0)
596         printf("PDG = %d (%s): generated: %d, reconstructed: %d, ratio: %f\n", 
597                (Int_t) fPIDParticles->GetBinCenter(i), 
598                pdgDB->GetParticle((Int_t) fPIDParticles->GetBinCenter(i))->GetName(), 
599                (Int_t) fPIDParticles->GetBinContent(i), 
600                (Int_t) fPIDTracks->GetBinContent(i), 
601                ((fPIDTracks->GetBinContent(i) > 0) ? fPIDParticles->GetBinContent(i) / fPIDTracks->GetBinContent(i) : -1));
602
603     delete pdgDB;
604     pdgDB = 0;
605   }
606
607   for (Int_t i=0; i<3; ++i)
608     if (fdNdEtaCorrectionVertexReco[i])
609       fdNdEtaCorrectionVertexReco[i]->Finish();
610
611   TFile* fout = TFile::Open("systematics.root", "RECREATE");
612
613   if (fEsdTrackCuts)
614     fEsdTrackCuts->SaveHistograms("esd_track_cuts");
615
616   if (fSecondaries)
617     fSecondaries->Write();
618
619   if (fSigmaVertex)
620     fSigmaVertex->Write();
621
622   for (Int_t i=0; i<4; ++i)
623     if (fdNdEtaCorrectionSpecies[i])
624       fdNdEtaCorrectionSpecies[i]->SaveHistograms();
625
626   for (Int_t i=0; i<3; ++i)
627     if (fdNdEtaCorrectionVertexReco[i])
628       fdNdEtaCorrectionVertexReco[i]->SaveHistograms();
629
630   fout->Write();
631   fout->Close();
632
633   if (fSecondaries)
634   {
635     new TCanvas;
636     fSecondaries->Draw("COLZ");
637   }
638
639   if (fSigmaVertex)
640   {
641     new TCanvas;
642     fSigmaVertex->Draw();
643   }
644 }