]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliFMDEventInspector.cxx
Change in centrality selection in AliFMDEventInspectror now it should correctly work...
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliFMDEventInspector.cxx
1 // 
2 // This class inspects the event 
3 //
4 // Input:
5 //   - AliESDFMD object possibly corrected for sharing
6 //
7 // Output:
8 //   - A histogram of v_z of events with triggers. 
9 //   - A histogram of v_z of events with vertex and triggers 
10 //   - A histogram of trigger counters 
11 // 
12 // Note, that these are added to the master output list 
13 //
14 // Corrections used: 
15 //   - None
16 //
17 #include "AliFMDEventInspector.h"
18 #include "AliLog.h"
19 #include "AliESDEvent.h"
20 #include "AliMultiplicity.h"
21 #include "AliAnalysisManager.h"
22 #include "AliMCEventHandler.h"
23 #include "AliInputEventHandler.h"
24 #include "AliTriggerAnalysis.h"
25 #include "AliPhysicsSelection.h"
26 #include "AliOADBPhysicsSelection.h"
27 #include "AliAODForwardMult.h"
28 #include "AliForwardUtil.h"
29 #include "AliCentrality.h"
30 #include <TH1.h>
31 #include <TList.h>
32 #include <TDirectory.h>
33 #include <TROOT.h>
34 #include <TParameter.h>
35 #include <iostream>
36 #include <iomanip>
37 #include "AliMCEvent.h"
38 #include "AliHeader.h"
39 #include "AliGenEventHeader.h"
40 #include "AliCollisionGeometry.h"
41 #include "AliVVZERO.h"
42
43 //====================================================================
44 const char* AliFMDEventInspector::fgkFolderName = "fmdEventInspector";
45
46 //____________________________________________________________________
47 AliFMDEventInspector::AliFMDEventInspector()
48   : TNamed(),
49     fHEventsTr(0), 
50     fHEventsTrVtx(0),
51     fHEventsAccepted(0),
52     fHEventsAcceptedXY(0),
53     fHTriggers(0),
54     fHTriggerCorr(0),
55     fHType(0),
56     fHWords(0),
57     fHCent(0),
58     fHCentVsQual(0),
59     fHStatus(0),
60     fLowFluxCut(1000),
61     fMaxVzErr(0.2),
62     fList(0),
63     fEnergy(0),
64     fField(999), 
65     fCollisionSystem(kUnknown),
66     fDebug(0),
67     fCentAxis(0),
68     fVtxAxis(10,-10,10),
69     fUseFirstPhysicsVertex(false),
70     fUseV0AND(false),
71     fMinPileupContrib(3), 
72   fMinPileupDistance(0.8),
73   fUseDisplacedVertices(false),
74   fDisplacedVertex(),
75   fCollWords(),
76   fBgWords(),
77   fCentMethod("V0M"),
78   fminCent(-1.0),
79   fmaxCent(-1.0),
80   fUsepA2012Vertex(false)                       
81 {
82   // 
83   // Constructor 
84   //
85   DGUARD(fDebug,1,"Default CTOR of AliFMDEventInspector");
86 }
87
88 //____________________________________________________________________
89 AliFMDEventInspector::AliFMDEventInspector(const char* name)
90   : TNamed(fgkFolderName, name),
91     fHEventsTr(0), 
92     fHEventsTrVtx(0), 
93     fHEventsAccepted(0),
94     fHEventsAcceptedXY(0),
95     fHTriggers(0),
96     fHTriggerCorr(0),
97     fHType(0),
98     fHWords(0),
99     fHCent(0),
100     fHCentVsQual(0),
101     fHStatus(0),
102     fLowFluxCut(1000),
103     fMaxVzErr(0.2),
104     fList(0),
105     fEnergy(0),
106     fField(999), 
107     fCollisionSystem(kUnknown),
108     fDebug(0),
109     fCentAxis(0),
110     fVtxAxis(10,-10,10),
111     fUseFirstPhysicsVertex(false),
112     fUseV0AND(false),
113     fMinPileupContrib(3), 
114     fMinPileupDistance(0.8),
115     fUseDisplacedVertices(false),
116   fDisplacedVertex(),
117   fCollWords(),
118   fBgWords(),
119   fCentMethod("V0M"),
120   fminCent(-1.0),
121   fmaxCent(-1.0),
122  fUsepA2012Vertex(false)                
123 {
124   // 
125   // Constructor 
126   // 
127   // Parameters:
128   //   name Name of object
129   //
130   DGUARD(fDebug,1,"Named CTOR of AliFMDEventInspector: %s", name);
131 }
132
133 //____________________________________________________________________
134 AliFMDEventInspector::AliFMDEventInspector(const AliFMDEventInspector& o)
135   : TNamed(o), 
136     fHEventsTr(o.fHEventsTr), 
137     fHEventsTrVtx(o.fHEventsTrVtx), 
138     fHEventsAccepted(o.fHEventsAccepted),
139     fHEventsAcceptedXY(o.fHEventsAcceptedXY),
140     fHTriggers(o.fHTriggers),
141     fHTriggerCorr(o.fHTriggerCorr),
142     fHType(o.fHType),
143     fHWords(o.fHWords),
144     fHCent(o.fHCent),
145     fHCentVsQual(o.fHCentVsQual),
146     fHStatus(o.fHStatus),
147     fLowFluxCut(o.fLowFluxCut),
148     fMaxVzErr(o.fMaxVzErr),
149     fList(o.fList),
150     fEnergy(o.fEnergy),
151     fField(o.fField), 
152     fCollisionSystem(o.fCollisionSystem),
153     fDebug(0),
154     fCentAxis(0),
155     fVtxAxis(o.fVtxAxis),
156     fUseFirstPhysicsVertex(o.fUseFirstPhysicsVertex),
157     fUseV0AND(o.fUseV0AND),
158     fMinPileupContrib(o.fMinPileupContrib), 
159     fMinPileupDistance(o.fMinPileupDistance),
160     fUseDisplacedVertices(o.fUseDisplacedVertices),
161     fDisplacedVertex(o.fDisplacedVertex),
162   fCollWords(),
163   fBgWords(),
164   fCentMethod(o.fCentMethod),
165   fminCent(o.fminCent),
166   fmaxCent(o.fmaxCent),
167   fUsepA2012Vertex(o.fUsepA2012Vertex)          
168 {
169   // 
170   // Copy constructor 
171   // 
172   // Parameters:
173   //   o Object to copy from 
174   //
175   DGUARD(fDebug,1,"Copy CTOR of AliFMDEventInspector");
176 }
177
178 //____________________________________________________________________
179 AliFMDEventInspector::~AliFMDEventInspector()
180 {
181   // 
182   // Destructor 
183   //
184   DGUARD(fDebug,1,"DTOR of AliFMDEventInspector");
185   // if (fList)         delete fList;
186 }
187 //____________________________________________________________________
188 AliFMDEventInspector&
189 AliFMDEventInspector::operator=(const AliFMDEventInspector& o)
190 {
191   // 
192   // Assignement operator
193   // 
194   // Parameters:
195   //   o Object to assign from 
196   // 
197   // Return:
198   //    Reference to this object
199   //
200   DGUARD(fDebug,3,"Assignment of AliFMDEventInspector");
201   if (&o == this) return *this; 
202   TNamed::operator=(o);
203   fHEventsTr         = o.fHEventsTr;
204   fHEventsTrVtx      = o.fHEventsTrVtx;
205   fHEventsAccepted   = o.fHEventsAccepted;
206   fHEventsAcceptedXY = o.fHEventsAcceptedXY;
207   fHTriggers         = o.fHTriggers;
208   fHTriggerCorr      = o.fHTriggerCorr;
209   fHType             = o.fHType;
210   fHWords            = o.fHWords;
211   fHCent             = o.fHCent;
212   fHCentVsQual       = o.fHCentVsQual;
213   fHStatus           = o.fHStatus;
214   fLowFluxCut        = o.fLowFluxCut;
215   fMaxVzErr          = o.fMaxVzErr;
216   fDebug             = o.fDebug;
217   fList              = (o.fList ? new TList : 0);
218   fEnergy            = o.fEnergy;
219   fField             = o.fField;
220   fCollisionSystem   = o.fCollisionSystem;
221   fVtxAxis.Set(o.fVtxAxis.GetNbins(), o.fVtxAxis.GetXmin(), 
222                o.fVtxAxis.GetXmax());
223   
224   fUseFirstPhysicsVertex = o.fUseFirstPhysicsVertex;
225   fUseV0AND              = o.fUseV0AND;
226   fMinPileupContrib      = o.fMinPileupContrib;
227   fMinPileupDistance     = o.fMinPileupDistance;
228   fUseDisplacedVertices  = o.fUseDisplacedVertices;
229   fDisplacedVertex       = o.fDisplacedVertex;
230   fCentMethod            = o.fCentMethod;
231   fminCent               = o.fminCent;
232   fmaxCent               = o.fmaxCent; 
233   fUsepA2012Vertex       =o.fUsepA2012Vertex;
234
235   if (fList) { 
236     fList->SetName(GetName());
237     if (fHEventsTr)    fList->Add(fHEventsTr);
238     if (fHEventsTrVtx) fList->Add(fHEventsTrVtx);
239     if (fHTriggers)    fList->Add(fHTriggers);
240     if (fHTriggerCorr) fList->Add(fHTriggerCorr);
241     if (fHType)        fList->Add(fHType);
242     if (fHWords)       fList->Add(fHWords);
243     if (fHCent)        fList->Add(fHCent);
244     if (fHCentVsQual)  fList->Add(fHCentVsQual);
245     if (fHStatus)      fList->Add(fHStatus);
246   }
247   return *this;
248 }
249
250 //____________________________________________________________________
251 void 
252 AliFMDEventInspector::SetCentralityMethod(ECentMethod m)
253 {
254   switch (m) { 
255   case kV0Multiplicity: fCentMethod = "VOM"; break; // VZERO multiplicity 
256   case kV0Amplitude:    fCentMethod = "V0A"; break; // VZERO amplitude    
257   case kV0Charge:       fCentMethod = "V0C"; break; // VZERO charge          
258   case kFMDRough:       fCentMethod = "FMD"; break; // FMD scaled energy l
259   case kNTracks:        fCentMethod = "TRK"; break; // Number of tracks   
260   case kLTracks:        fCentMethod = "TKL"; break; // Number of tracks   
261   case kCL0:            fCentMethod = "CL0"; break; //               
262   case kCL1:            fCentMethod = "CL1"; break; //               
263   case kCND:            fCentMethod = "CND"; break; //               
264   case kNParticles:     fCentMethod = "NPA"; break; // Neutral particles  
265   case kNeutrons:       fCentMethod = "ZNA"; break; // ZDC neutron amplitu
266   case kV0vsFMD:        fCentMethod = "V0MvsFMD"; break; // VZERO versus FMD   
267   case kV0vsNTracks:    fCentMethod = "TKLvsVOM"; break; // Tracks versus VZERO
268   case kZEMvsZDC:       fCentMethod = "ZEMvsZDC"; break; // ZDC              
269   default:              fCentMethod = "V0M"; break;
270   }
271 }
272
273 //____________________________________________________________________
274 Bool_t
275 AliFMDEventInspector::FetchHistograms(const TList* d, 
276                                       TH1I*& hEventsTr, 
277                                       TH1I*& hEventsTrVtx, 
278                                       TH1I*& hEventsAcc,
279                                       TH1I*& hTriggers) const
280 {
281   // 
282   // Fetch our histograms from the passed list 
283   // 
284   // Parameters:
285   //   d             Input
286   //   hEventsTr     On return, pointer to histogram, or null
287   //   hEventsTrVtx  On return, pointer to histogram, or null
288   //   hTriggers     On return, pointer to histogram, or null
289   // 
290   // Return:
291   //    true on success, false otherwise 
292   //
293   DGUARD(fDebug,3,"Fetch histograms in AliFMDEventInspector");
294   hEventsTr    = 0;
295   hEventsTrVtx = 0;
296   hEventsAcc   = 0;
297   hTriggers    = 0;  
298   TList* dd    = dynamic_cast<TList*>(d->FindObject(GetName()));
299   if (!dd) return kFALSE;
300   
301   hEventsTr    = dynamic_cast<TH1I*>(dd->FindObject("nEventsTr"));
302   hEventsTrVtx = dynamic_cast<TH1I*>(dd->FindObject("nEventsTrVtx"));
303   hEventsAcc   = dynamic_cast<TH1I*>(dd->FindObject("nEventsAccepted"));
304   hTriggers    = dynamic_cast<TH1I*>(dd->FindObject("triggers"));
305
306   if (!hEventsTr    || 
307       !hEventsTrVtx || 
308       !hEventsAcc   ||
309       !hTriggers) return kFALSE;
310   return kTRUE;
311 }
312 //____________________________________________________________________
313 void
314 AliFMDEventInspector::CacheConfiguredTriggerClasses(TList& cache, 
315                                                     const TList* classes,
316                                                     AliOADBPhysicsSelection* o)
317 {
318   TIter nextClass(classes);
319   TObjString* trigClass = 0;
320   // Loop over all trigger classes.  Trigger classes have the format 
321   //
322   //   class          := positive_words SPACE(s) negative_words 
323   //   positive_words := 
324   //                  |  '+' words
325   //   negative_words := 
326   //                  |  '-' words
327   //   words          := word 
328   //                  |  word ',' words 
329   //   
330   while ((trigClass = static_cast<TObjString*>(nextClass()))) {
331     // Tokenize on space to get positive and negative parts 
332     TString     side   = o->GetBeamSide(trigClass->String());
333     TObjArray*  parts  = trigClass->String().Tokenize(" ");
334     TObjString* part   = 0;
335     TIter       nextPart(parts);
336     while ((part = static_cast<TObjString*>(nextPart()))) {
337       // We only care about the positive ones 
338       if (part->GetName()[0] != '+') continue;
339       part->String().Remove(0,1);
340         
341       // Tokenize on a comma to get the words 
342       TObjArray*  words = part->String().Tokenize(",");
343       TObjString* word  = 0;
344       TIter       nextWord(words);
345       while ((word = static_cast<TObjString*>(nextWord()))) {
346         TNamed* store = new TNamed(word->String(), side);
347         cache.Add(store);
348         DMSG(fDebug,3,"Caching %s trigger word %s", 
349              store->GetTitle(), store->GetName());
350       } // while (word)
351       delete words;
352     }
353     delete parts;
354   }
355 }
356
357 //____________________________________________________________________
358 void
359 AliFMDEventInspector::SetupForData(const TAxis& vtxAxis)
360 {
361   // 
362   // Initialize the object - this is called on the first seen event. 
363   // 
364   // Parameters:
365   //   vtxAxis Vertex axis in use 
366   //
367   DGUARD(fDebug,1,"Initialize in AliFMDEventInspector");
368
369   AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
370
371   // Get the input handler - should always be there 
372   AliInputEventHandler* ih = 
373     static_cast<AliInputEventHandler*>(am->GetInputEventHandler());
374   if (!ih) { 
375     AliWarning("No input handler");
376     return;
377   }
378   // Get the physics selection - should always be there 
379   AliPhysicsSelection* ps = 
380     static_cast<AliPhysicsSelection*>(ih->GetEventSelection());
381   if (!ps) {
382     AliWarning("No physics selection");
383     return;
384   }
385   // Get the configured triggers
386   AliOADBPhysicsSelection* oadb = 
387     const_cast<AliOADBPhysicsSelection*>(ps->GetOADBPhysicsSelection());
388   if (!oadb) {
389     AliWarning("No OADB physics selection object");
390     return;
391   }
392   // Get the configured trigger words from the physics selection 
393   const TList* collTriggClasses = ps->GetCollisionTriggerClasses();
394   const TList* bgTriggClasses   = ps->GetBGTriggerClasses();
395   if (!collTriggClasses) { 
396     AliWarning("No configured collision trigger classes");
397     return;
398   }
399   if (!bgTriggClasses) { 
400     AliWarning("No configured background trigger classes");
401     return;
402   }
403   CacheConfiguredTriggerClasses(fCollWords, collTriggClasses, oadb);
404   CacheConfiguredTriggerClasses(fBgWords,   bgTriggClasses,   oadb);
405   // fCollWords.ls();
406   // fBgWords.ls();
407   
408   
409   // -1.5 -0.5 0.5 1.5 ... 89.5 ... 100.5
410   // ----- 92 number --------- ---- 1 ---
411   TArrayD limits(93);
412   for (Int_t i = 0; i < 92; i++) limits[i] = -1.5 + i;
413   limits[92] = 100.5;
414
415   fVtxAxis.Set(vtxAxis.GetNbins(), vtxAxis.GetXmin(), vtxAxis.GetXmax());
416   
417   fCentAxis  = new TAxis(limits.GetSize()-1, limits.GetArray());
418   fHEventsTr = new TH1I("nEventsTr", "Number of events w/trigger", 
419                         4*vtxAxis.GetNbins(), 
420                         2*vtxAxis.GetXmin(), 
421                         2*vtxAxis.GetXmax());
422   fHEventsTr->SetXTitle("v_{z} [cm]");
423   fHEventsTr->SetYTitle("# of events");
424   fHEventsTr->SetFillColor(kRed+1);
425   fHEventsTr->SetFillStyle(3001);
426   fHEventsTr->SetDirectory(0);
427   // fHEventsTr->Sumw2();
428   fList->Add(fHEventsTr);
429
430   fHEventsTrVtx = static_cast<TH1I*>(fHEventsTr->Clone("nEventsTrVtx"));
431   fHEventsTrVtx->SetTitle("Number of events w/trigger and vertex"); 
432   fHEventsTrVtx->SetFillColor(kBlue+1);
433   fHEventsTrVtx->SetDirectory(0);
434   // fHEventsTrVtx->Sumw2();
435   fList->Add(fHEventsTrVtx);
436
437   fHEventsAccepted = new TH1I("nEventsAccepted", 
438                               "Number of events  w/trigger and vertex in range",
439                               2*vtxAxis.GetNbins(), 
440                               2*vtxAxis.GetXmin(), 
441                               2*vtxAxis.GetXmax());
442   fHEventsAccepted->SetXTitle("v_{z} [cm]");
443   fHEventsAccepted->SetYTitle("# of events");
444   fHEventsAccepted->SetFillColor(kGreen+1);
445   fHEventsAccepted->SetFillStyle(3001);
446   fHEventsAccepted->SetDirectory(0);
447   // fHEventsAccepted->Sumw2();
448   fList->Add(fHEventsAccepted);
449                               
450   fHEventsAcceptedXY = new TH2D("nEventsAcceptedXY", 
451                               "XY vertex w/trigger and Z vertex in range",
452                                 1000,-1,1,1000,-1,1);
453   
454   fHEventsAcceptedXY->SetXTitle("v_{x} [cm]");
455   fHEventsAcceptedXY->SetYTitle("v_{y} [cm]");
456   fHEventsAcceptedXY->SetDirectory(0);
457   // fHEventsAccepted->Sumw2();
458   fList->Add(fHEventsAcceptedXY);
459
460   
461   fHTriggers = new TH1I("triggers", "Triggers", kOffline+1, 0, kOffline+1);
462   fHTriggers->SetFillColor(kRed+1);
463   fHTriggers->SetFillStyle(3001);
464   fHTriggers->SetStats(0);
465   fHTriggers->SetDirectory(0);
466
467   fHTriggerCorr = new TH2I("triggerCorr", "Trigger correlation", 
468                            kOffline+1, 0, kOffline+1, 
469                            kOffline+1, 0, kOffline+1);
470   fHTriggerCorr->SetStats(0);
471   fHTriggerCorr->SetDirectory(0);
472   fHTriggerCorr->SetXTitle("Requirement");
473   fHTriggerCorr->SetYTitle("Companion");
474
475   Int_t binNum[] = { kInel   +1,
476                      kInelGt0+1,
477                      kNSD    +1,
478                      kV0AND  +1,
479                      kEmpty  +1,
480                      kA      +1,
481                      kB      +1,
482                      kC      +1,
483                      kE      +1,
484                      kPileUp +1,
485                      kMCNSD  +1,
486                      kOffline+1 };
487   const char* binLbl[] = { "INEL",       
488                            "INEL>0",
489                            "NSD",        
490                            "VOAND",
491                            "Empty",      
492                            "A",  
493                            "B",  
494                            "C",  
495                            "E",  
496                            "Pileup",
497                            "NSD_{MC}", 
498                            "Offline" };
499   for (Int_t i = 0; i < kOffline+1; i++) {
500     fHTriggers->GetXaxis()->SetBinLabel(binNum[i], binLbl[i]);
501     fHTriggerCorr->GetXaxis()->SetBinLabel(binNum[i], binLbl[i]);
502     fHTriggerCorr->GetYaxis()->SetBinLabel(binNum[i], binLbl[i]);
503   }
504   fList->Add(fHTriggers);
505   fList->Add(fHTriggerCorr);
506   
507
508   fHType = new TH1I("type", Form("Event type (cut: SPD mult>%d)", 
509                                  fLowFluxCut), 2, -.5, 1.5);
510   fHType->SetFillColor(kRed+1);
511   fHType->SetFillStyle(3001);
512   fHType->SetStats(0);
513   fHType->SetDirectory(0);
514   fHType->GetXaxis()->SetBinLabel(1,"Low-flux");
515   fHType->GetXaxis()->SetBinLabel(2,"High-flux");
516   fList->Add(fHType);
517
518 #if 0 
519   // This histogram disabled as it causes problems in the merge 
520   fHWords = new TH1I("words", "Trigger words seen", 1, 0, 0); 
521   fHWords->SetFillColor(kBlue+1);
522   fHWords->SetFillStyle(3001);
523   fHWords->SetStats(0);
524   fHWords->SetDirectory(0);
525   fHWords->SetBit(TH1::kCanRebin);
526   fList->Add(fHWords);
527 #endif
528
529   fHCent = new TH1F("cent", "Centrality", limits.GetSize()-1,limits.GetArray());
530   fHCent->SetFillColor(kBlue+1);
531   fHCent->SetFillStyle(3001);
532   fHCent->SetStats(0);
533   fHCent->SetDirectory(0);
534   fHCent->SetXTitle("Centrality [%]");
535   fHCent->SetYTitle("Events");
536   fList->Add(fHCent);
537
538   fHCentVsQual = new TH2F("centVsQuality", "Quality vs Centrality", 
539                           5, 0, 5, limits.GetSize()-1, limits.GetArray());
540   fHCentVsQual->SetXTitle("Quality");
541   fHCentVsQual->SetYTitle("Centrality [%]");
542   fHCentVsQual->SetZTitle("Events");
543   fHCentVsQual->GetXaxis()->SetBinLabel(1, "OK");
544   fHCentVsQual->GetXaxis()->SetBinLabel(2, "Outside v_{z} cut");
545   fHCentVsQual->GetXaxis()->SetBinLabel(3, "V0 vs SPD outlier");
546   fHCentVsQual->GetXaxis()->SetBinLabel(4, "V0 vs TPC outlier");
547   fHCentVsQual->GetXaxis()->SetBinLabel(5, "V0 vs ZDC outlier");
548   fHCentVsQual->SetDirectory(0);
549   fList->Add(fHCentVsQual);
550
551   fHStatus = new TH1I("status", "Status", 7, 1, 8);
552   fHStatus->SetFillColor(kRed+1);
553   fHStatus->SetFillStyle(3001);
554   fHStatus->SetStats(0);
555   fHStatus->SetDirectory(0);
556   fHStatus->GetXaxis()->SetBinLabel(1, "OK");
557   fHStatus->GetXaxis()->SetBinLabel(2, "No event");
558   fHStatus->GetXaxis()->SetBinLabel(3, "No triggers");
559   fHStatus->GetXaxis()->SetBinLabel(4, "No SPD");
560   fHStatus->GetXaxis()->SetBinLabel(5, "No FMD");
561   fHStatus->GetXaxis()->SetBinLabel(6, "No vertex");
562   fHStatus->GetXaxis()->SetBinLabel(7, "Bad vertex");
563   fList->Add(fHStatus);
564 }
565
566 //____________________________________________________________________
567 void
568 AliFMDEventInspector::StoreInformation(Int_t runNo)
569 {
570   // Write TNamed objects to output list containing information about
571   // the running conditions 
572   DGUARD(fDebug,2,"Store information from AliFMDEventInspector");
573   if (!fList) return;
574
575   
576   fList->Add(AliForwardUtil::MakeParameter("sys", fCollisionSystem));
577   fList->Add(AliForwardUtil::MakeParameter("sNN", fEnergy));
578   fList->Add(AliForwardUtil::MakeParameter("field", fField));
579   fList->Add(AliForwardUtil::MakeParameter("runNo", runNo));
580   fList->Add(AliForwardUtil::MakeParameter("lowFlux", fLowFluxCut));
581   fList->Add(AliForwardUtil::MakeParameter("fpVtx",fUseFirstPhysicsVertex));
582   fList->Add(AliForwardUtil::MakeParameter("v0and",fUseV0AND));
583   fList->Add(AliForwardUtil::MakeParameter("nPileUp", fMinPileupContrib));
584   fList->Add(AliForwardUtil::MakeParameter("dPileup", fMinPileupDistance));
585   fList->Add(AliForwardUtil::MakeParameter("alirootRev", 
586                                            AliForwardUtil::AliROOTRevision()));
587   fList->Add(AliForwardUtil::MakeParameter("alirootBranch", 
588                                            AliForwardUtil::AliROOTBranch()));
589
590 }
591
592 //____________________________________________________________________
593 void
594 AliFMDEventInspector::CreateOutputObjects(TList* dir)
595 {
596   // 
597   // Define the output histograms.  These are put in a sub list of the
598   // passed list.   The histograms are merged before the parent task calls 
599   // AliAnalysisTaskSE::Terminate 
600   // 
601   //   dir Directory to add to 
602   //
603   DGUARD(fDebug,1,"Define output from AliFMDEventInspector");
604   fList = new TList;
605   fList->SetName(GetName());
606   fList->SetOwner();
607   dir->Add(fList);
608 }
609
610 //____________________________________________________________________
611 UInt_t
612 AliFMDEventInspector::Process(const AliESDEvent* event, 
613                               UInt_t&            triggers,
614                               Bool_t&            lowFlux,
615                               UShort_t&          ivz, 
616                               TVector3&          ip,
617                               Double_t&          cent,
618                               UShort_t&          nClusters)
619 {
620   // 
621   // Process the event 
622   // 
623   // Parameters:
624   //   event     Input event 
625   //   triggers  On return, the triggers fired 
626   //   lowFlux   On return, true if the event is considered a low-flux 
627   //                  event (according to the setting of fLowFluxCut) 
628   //   ivz       On return, the found vertex bin (1-based).  A zero
629   //                  means outside of the defined vertex range
630   //   vz        On return, the z position of the interaction
631   //   cent      On return, the centrality - if not available < 0
632   // 
633   // Return:
634   //    0 (or kOk) on success, otherwise a bit mask of error codes 
635   //
636   DGUARD(fDebug,1,"Process event in AliFMDEventInspector"); 
637   // --- Check that we have an event ---------------------------------
638   if (!event) { 
639     AliWarning("No ESD event found for input event");
640     fHStatus->Fill(2);
641     return kNoEvent;
642   }
643
644   // --- Read trigger information from the ESD and store in AOD object
645   if (!ReadTriggers(*event, triggers, nClusters)) { 
646     if (fDebug > 2) {
647       AliWarning("Failed to read triggers from ESD"); }
648     fHStatus->Fill(3);
649     return kNoTriggers;
650   }
651
652   // --- Check if this is a high-flux event --------------------------
653   const AliMultiplicity* testmult = event->GetMultiplicity();
654   if (!testmult) {
655     if (fDebug > 3) {
656       AliWarning("No central multiplicity object found"); }
657   }
658   else 
659     lowFlux = testmult->GetNumberOfTracklets() < fLowFluxCut;
660
661   fHType->Fill(lowFlux ? 0 : 1);
662
663   // --- Process satellite event information is requested ------------
664   if (fUseDisplacedVertices) { 
665     if (!fDisplacedVertex.Process(event)) 
666       AliWarning("Failed to process satellite event");
667   }
668   
669  
670   // --- Get the interaction point -----------------------------------
671   Bool_t vzOk = ReadVertex(*event, ip);
672   fHEventsTr->Fill(ip.Z());
673   if (!vzOk) { 
674     if (fDebug > 3) {
675       AliWarning("Failed to read vertex from ESD"); }
676     fHStatus->Fill(6);
677     return kNoVertex;
678   }
679
680    // --- Read centrality information 
681   cent          = -10;
682   UShort_t qual = 0;
683   if (!ReadCentrality(*event, cent, qual)) {
684     if (fDebug > 3) 
685       AliWarning("Failed to get centrality");
686   }
687   // --- check centrality cut
688  
689   if(fminCent>-0.0001&&cent<fminCent)
690         return  kNoEvent; 
691   if(fmaxCent>-0.0001&&cent>fmaxCent)
692          return  kNoEvent; 
693   fHCent->Fill(cent);
694   if (qual == 0) fHCentVsQual->Fill(0., cent);
695   else { 
696     for (UShort_t i = 0; i < 4; i++) 
697       if (qual & (1 << i)) fHCentVsQual->Fill(Double_t(i+1), cent);
698   }
699         
700
701   fHEventsTrVtx->Fill(ip.Z());
702   
703   // --- Get the vertex bin ------------------------------------------
704   ivz = fVtxAxis.FindBin(ip.Z());
705   if (ivz <= 0 || ivz > fVtxAxis.GetNbins()) { 
706     if (fDebug > 3) {
707       AliWarning(Form("Vertex @ %f outside of range [%f,%f]", 
708                       ip.Z(), fVtxAxis.GetXmin(), fVtxAxis.GetXmax())); 
709     }
710     ivz = 0;
711     fHStatus->Fill(7);
712     return kBadVertex;
713   }
714   fHEventsAccepted->Fill(ip.Z());
715   fHEventsAcceptedXY->Fill(ip.X(),ip.Y());
716   
717   // --- Check the FMD ESD data --------------------------------------
718   if (!event->GetFMDData()) { 
719     if (fDebug > 3) {
720       AliWarning("No FMD data found in ESD"); }
721     fHStatus->Fill(5);
722     return kNoFMD;
723   }
724
725   fHStatus->Fill(1);
726   return kOk;
727 }
728
729 //____________________________________________________________________
730 Bool_t
731 AliFMDEventInspector::ReadCentrality(const AliESDEvent& esd, 
732                                      Double_t& cent, 
733                                      UShort_t& qual) const
734 {
735   // 
736   // Read centrality from event 
737   // 
738   // Parameters:
739   //    esd  Event 
740   //    cent On return, the centrality or negative if not found
741   // 
742   // Return:
743   //    False on error, true otherwise 
744   //
745   DGUARD(fDebug,2,"Read the centrality in AliFMDEventInspector");
746
747   if(fUseDisplacedVertices) {
748     Double_t zvtx = fDisplacedVertex.GetVertexZ();
749     qual          = 1;
750     if(TMath::Abs(zvtx) < 999) {
751       cent = fDisplacedVertex.GetCentralityPercentile();
752       qual = 0;
753     }
754     return true;
755   }
756   
757   cent = -1;
758   qual = 0;
759   AliCentrality* centObj = const_cast<AliESDEvent&>(esd).GetCentrality();
760   if (!centObj)  return true;
761
762   cent = centObj->GetCentralityPercentile(fCentMethod);  
763   qual = centObj->GetQuality();
764
765   return true;
766 }
767
768 //____________________________________________________________________
769 Bool_t
770 AliFMDEventInspector::CheckpAExtraV0(const AliESDEvent& esd) const
771 {
772   if (fCollisionSystem != AliForwardUtil::kPPb) return true;
773
774    AliVVZERO* esdV0 = esd.GetVZEROData();
775    if ((esdV0->GetV0ADecision()!=1) || (esdV0->GetV0CDecision()!=1)) 
776      return false;
777    return true;
778 }
779
780 //____________________________________________________________________
781 Bool_t
782 AliFMDEventInspector::ReadTriggers(const AliESDEvent& esd, UInt_t& triggers,
783                                    UShort_t& nClusters)
784 {
785   // 
786   // Read the trigger information from the ESD event 
787   // 
788   // Parameters:
789   //   esd        ESD event 
790   //   triggers   On return, contains the trigger bits 
791   // 
792   // Return:
793   //    @c true on success, @c false otherwise 
794   //
795   DGUARD(fDebug,2,"Read the triggers in AliFMDEventInspector");
796   triggers = 0;
797
798   // Get the analysis manager - should always be there 
799   AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
800   DMSG(fDebug,10,"Got analysis manager %p", am);
801   if (!am) { 
802     AliWarning("No analysis manager defined!");
803     return kFALSE;
804   }
805
806   // Get the input handler - should always be there 
807   AliInputEventHandler* ih = 
808     static_cast<AliInputEventHandler*>(am->GetInputEventHandler());
809   DMSG(fDebug,10,"Got input handler %p", ih);
810   if (!ih) { 
811     AliWarning("No input handler");
812     return kFALSE;
813   }
814
815   // Check if this is a collision candidate (MB)
816   // Note, that we should use the value cached in the input 
817   // handler rather than calling IsCollisionCandiate directly 
818   // on the AliPhysicsSelection obejct.  If we called the latter
819   // then the AliPhysicsSelection object would overcount by a 
820   // factor of 2! :-(
821   Bool_t  offline  = ih->IsEventSelected();
822   Bool_t  fastonly = (ih->IsEventSelected() & AliVEvent::kFastOnly);
823   TString trigStr  = esd.GetFiredTriggerClasses();
824
825   if (fHWords) fHWords->Fill(trigStr.Data(), 1);
826   
827   if(fUseDisplacedVertices) {
828     DMSG(fDebug,3,"Using displaced vertex stuff");
829     if (TMath::Abs(fDisplacedVertex.GetVertexZ()) >= 999) offline = false;
830   }
831   
832   if (CheckFastPartition(fastonly))     offline = false;
833   if (offline && CheckCosmics(trigStr)) offline = false;
834   // if (!CheckpAExtraV0(esd))             offline = false;
835
836   DMSG(fDebug,2,"Event is %striggered by off-line", offline ? "" : "NOT ");
837
838   if (offline) {
839     triggers |= AliAODForwardMult::kOffline;
840     triggers |= AliAODForwardMult::kInel;
841     if (!fHTriggers) { 
842       AliWarning("Histogram of triggers not defined - has init been called");
843       return false;
844     }
845     // fHTriggers->Fill(kOffline+0.5);
846     
847     CheckINELGT0(esd, nClusters, triggers);
848   }
849   
850   CheckNSD(esd,triggers);
851   CheckPileup(esd, triggers);
852   CheckEmpty(trigStr, triggers);
853   // if (CheckPileup(esd, triggers)) fHTriggers->Fill(kPileUp+.5);
854   // if (CheckEmpty(trigStr, triggers)) fHTriggers->Fill(kEmpty+.5);
855
856   CheckWords(esd, triggers);
857
858 #if 0
859   // Now check - if we have a collision - for offline triggers and
860   // fill histogram.
861   if (triggers & AliAODForwardMult::kB) {
862     fHTriggers->Fill(kB+.5);
863     if (triggers & AliAODForwardMult::kInel) 
864       fHTriggers->Fill(kInel+.5);
865     
866     if (triggers & AliAODForwardMult::kInelGt0)
867       fHTriggers->Fill(kInelGt0+.5);
868     
869     if (triggers & AliAODForwardMult::kNSD)
870       fHTriggers->Fill(kNSD+.5);
871
872     if (triggers & AliAODForwardMult::kV0AND)
873       fHTriggers->Fill(kV0AND+.5);
874   }
875   if (triggers & AliAODForwardMult::kA) fHTriggers->Fill(kA+.5);
876   if (triggers & AliAODForwardMult::kC) fHTriggers->Fill(kC+.5);
877   if (triggers & AliAODForwardMult::kE) fHTriggers->Fill(kE+.5);
878 #endif
879 #define TEST_TRIG_BIN(RET,BIN,TRIGGERS) \
880   do { switch (BIN) { \
881     case kInel:    RET = triggers & AliAODForwardMult::kInel;    break; \
882     case kInelGt0: RET = triggers & AliAODForwardMult::kInelGt0; break; \
883     case kNSD:     RET = triggers & AliAODForwardMult::kNSD;     break; \
884     case kV0AND:   RET = triggers & AliAODForwardMult::kV0AND;   break; \
885     case kEmpty:   RET = triggers & AliAODForwardMult::kEmpty;   break; \
886     case kA:       RET = triggers & AliAODForwardMult::kA;       break; \
887     case kB:       RET = triggers & AliAODForwardMult::kB;       break; \
888     case kC:       RET = triggers & AliAODForwardMult::kC;       break; \
889     case kE:       RET = triggers & AliAODForwardMult::kE;       break; \
890     case kPileUp:  RET = triggers & AliAODForwardMult::kPileUp;  break; \
891     case kMCNSD:   RET = triggers & AliAODForwardMult::kMCNSD;   break; \
892     case kOffline: RET = triggers & AliAODForwardMult::kOffline; break; \
893     default:       RET = false; } } while(false)
894       
895   
896   for (Int_t i = 0; i < kOffline+1; i++) { 
897     Bool_t hasX = false;
898     TEST_TRIG_BIN(hasX, i, triggers);
899     if (!hasX) continue;
900     fHTriggers->Fill(i+.5);
901     for (Int_t j = 0; j < kOffline+1; j++) { 
902       Bool_t hasY = false;
903       TEST_TRIG_BIN(hasY, j, triggers);
904       if (!hasY) continue;
905       
906       fHTriggerCorr->Fill(i+.5, j+.5);
907     }
908   }
909   return kTRUE;
910 }
911
912 //____________________________________________________________________
913 Bool_t
914 AliFMDEventInspector::CheckFastPartition(bool fastonly) const
915 {
916   // For the 2.76 TeV p+p run, the FMD ran in the slow partition 
917   // so it received no triggers from the fast partition. Therefore
918   // the fast triggers are removed here but not for MC where all 
919   // triggers are fast.
920   if (TMath::Abs(fEnergy - 2750.) > 20) return false;
921   if (fCollisionSystem != AliForwardUtil::kPP) return false;
922   if (fastonly)
923     DMSG(fDebug,1,"Fast trigger in pp @ sqrt(s)=2.76TeV removed");
924
925   return fastonly;
926 }
927
928 //____________________________________________________________________
929 Bool_t
930 AliFMDEventInspector::CheckCosmics(const TString& trigStr) const
931 {
932   // MUON triggers are not strictly minimum bias (MB) so they are
933   // removed (HHD)
934   if(trigStr.Contains("CMUS1")) {
935     DMSG(fDebug,1,"Cosmic trigger ins't min-bias, removed");
936     return true;
937   }
938   return false;
939 }
940
941 //____________________________________________________________________
942 Bool_t
943 AliFMDEventInspector::CheckINELGT0(const AliESDEvent& esd, 
944                                    UShort_t& nClusters, 
945                                    UInt_t& triggers) const
946 {
947   nClusters = 0;
948
949   // If this is inel, see if we have a tracklet 
950   const AliMultiplicity* spdmult = esd.GetMultiplicity();
951   if (!spdmult) {
952     AliWarning("No SPD multiplicity");
953     return false;
954   }
955
956   // Check if we have one or more tracklets 
957   // in the range -1 < eta < 1 to set the INEL>0 
958   // trigger flag. 
959   // 
960   // Also count tracklets as a single cluster 
961   Int_t n = spdmult->GetNumberOfTracklets();
962   for (Int_t j = 0; j < n; j++) { 
963     if(TMath::Abs(spdmult->GetEta(j)) < 1) { 
964       triggers |= AliAODForwardMult::kInelGt0;
965       nClusters++;
966     }
967   }
968   n = spdmult->GetNumberOfSingleClusters();
969   for (Int_t j = 0; j < n; j++) { 
970     Double_t eta = -TMath::Log(TMath::Tan(spdmult->GetThetaSingle(j)/2.));
971     if (TMath::Abs(eta) < 1) nClusters++;
972   }
973   if (nClusters > 0) triggers |= AliAODForwardMult::kNClusterGt0;
974
975   return triggers & AliAODForwardMult::kNClusterGt0;
976 }
977
978 //____________________________________________________________________
979 Bool_t
980 AliFMDEventInspector::CheckNSD(const AliESDEvent& esd, UInt_t& triggers) const
981 {
982   // Analyse some trigger stuff 
983   AliTriggerAnalysis ta;
984   if (ta.IsOfflineTriggerFired(&esd, AliTriggerAnalysis::kV0AND)) {
985     triggers |= AliAODForwardMult::kV0AND;
986     if (fUseV0AND) 
987       triggers |= AliAODForwardMult::kNSD;
988   }
989   if (ta.IsOfflineTriggerFired(&esd, AliTriggerAnalysis::kNSD1)) 
990     triggers |= AliAODForwardMult::kNSD;
991   return triggers & AliAODForwardMult::kNSD;
992 }
993 //____________________________________________________________________
994 Bool_t
995 AliFMDEventInspector::CheckPileup(const AliESDEvent& esd, 
996                                   UInt_t& triggers) const
997 {
998   // Check for multiple vertices (pile-up) with at least 3
999   // contributors and at least 0.8cm from the primary vertex
1000   if(fCollisionSystem != AliForwardUtil::kPP) return false;
1001
1002   Bool_t pileup =  esd.IsPileupFromSPD(fMinPileupContrib,fMinPileupDistance);
1003   if (pileup) triggers |= AliAODForwardMult::kPileUp;
1004   return pileup;
1005 }
1006
1007 //____________________________________________________________________
1008 Bool_t
1009 AliFMDEventInspector::CheckEmpty(const TString& trigStr, UInt_t& triggers) const
1010 {
1011   if (trigStr.Contains("CBEAMB-ABCE-NOPF-ALL")) {
1012     triggers |= AliAODForwardMult::kEmpty;
1013     return true;
1014   }
1015   return false;
1016 }
1017 //____________________________________________________________________
1018 Bool_t
1019 AliFMDEventInspector::CheckWords(const AliESDEvent& esd, UInt_t& triggers) const
1020 {
1021   TObject* word = 0;
1022   TIter nextColl(&fCollWords);
1023   while ((word = nextColl())) {
1024     DMSG(fDebug,10,"Checking if %s trigger %s is fired", 
1025          word->GetTitle(), word->GetName());
1026     if (!esd.IsTriggerClassFired(word->GetName())) continue;
1027
1028     TString beamSide = word->GetTitle();
1029     DMSG(fDebug,10,"Found it - this is a %s trigger", beamSide.Data());
1030
1031     if (!beamSide.EqualTo("B")) continue;
1032     triggers |= AliAODForwardMult::kB;
1033     break; // No more to do here
1034   }
1035   TIter nextBg(&fBgWords);
1036   UInt_t all = (AliAODForwardMult::kA | 
1037                 AliAODForwardMult::kC | 
1038                 AliAODForwardMult::kE);
1039   while ((word = nextBg())) {
1040     DMSG(fDebug,10,"Checking if %s trigger %s is fired", 
1041          word->GetTitle(), word->GetName());
1042     if (!esd.IsTriggerClassFired(word->GetName())) continue;
1043     
1044     TString beamSide = word->GetTitle();
1045     DMSG(fDebug,10,"Found it - this is a %s trigger", beamSide.Data());
1046
1047     if (beamSide.Contains("A")) triggers |= AliAODForwardMult::kA;
1048     if (beamSide.Contains("C")) triggers |= AliAODForwardMult::kC;
1049     if (beamSide.Contains("E")) triggers |= AliAODForwardMult::kE;
1050
1051     if ((triggers & all) == all) break; // No more to do
1052   }
1053   return true;
1054 }
1055
1056
1057 //____________________________________________________________________
1058 Bool_t
1059 AliFMDEventInspector::ReadVertex(const AliESDEvent& esd, TVector3& ip)
1060 {
1061   // 
1062   // Read the vertex information from the ESD event 
1063   // 
1064   // Parameters:
1065   //   esd  ESD event 
1066   //   vz   On return, the vertex Z position 
1067   // 
1068   // Return:
1069   //    @c true on success, @c false otherwise 
1070   //
1071   DGUARD(fDebug,2,"Read the vertex in AliFMDEventInspector");
1072   ip.SetXYZ(1024, 1024, 0);
1073   
1074   if(fUseDisplacedVertices) {
1075     Double_t zvtx = fDisplacedVertex.GetVertexZ();
1076       
1077     if(TMath::Abs(zvtx) < 999) {
1078       ip.SetZ(zvtx);
1079       return true;
1080     }
1081     return false;
1082   }
1083
1084   if(fUseFirstPhysicsVertex) return CheckPWGUDVertex(esd, ip);
1085   
1086   if(fUsepA2012Vertex) return CheckpA2012Vertex(esd,ip);        
1087   
1088   return CheckVertex(esd, ip);
1089 }
1090
1091 //____________________________________________________________________
1092 Bool_t
1093 AliFMDEventInspector::CheckPWGUDVertex(const AliESDEvent& esd, 
1094                                        TVector3& ip)  const
1095 {
1096   // This is the code used by the 1st physics people 
1097   const AliESDVertex* vertex    = esd.GetPrimaryVertex();
1098   if (!vertex  || !vertex->GetStatus()) {
1099     DMSG(fDebug,2,"No primary vertex (%p) or bad status %d", 
1100          vertex, (vertex ? vertex->GetStatus() : -1));
1101     return false;
1102   }
1103   const AliESDVertex* vertexSPD = esd.GetPrimaryVertexSPD();
1104   if (!vertexSPD || !vertexSPD->GetStatus()) {
1105     DMSG(fDebug,2,"No primary SPD vertex (%p) or bad status %d", 
1106          vertexSPD, (vertexSPD ? vertexSPD->GetStatus() : -1));
1107     return false;
1108   }
1109     
1110   // if vertex is from SPD vertexZ, require more stringent cuts 
1111   if (vertex->IsFromVertexerZ()) { 
1112     if (vertex->GetDispersion() > fMaxVzErr || 
1113         vertex->GetZRes() > 1.25 * fMaxVzErr) {
1114       DMSG(fDebug,2,"Dispersion %f > %f or resolution %f > %f",
1115            vertex->GetDispersion(), fMaxVzErr,
1116            vertex->GetZRes(), 1.25 * fMaxVzErr);
1117       return false;
1118     }
1119   }
1120   ip.SetZ(vertex->GetZ());
1121   
1122   if(!vertex->IsFromVertexerZ()) {
1123     ip.SetX(vertex->GetX());
1124     ip.SetY(vertex->GetY());
1125   }
1126   return true;
1127 }
1128 //
1129 Bool_t AliFMDEventInspector::CheckpA2012Vertex(const AliESDEvent& esd, 
1130                                        TVector3& ip)  const
1131 {      
1132       const AliESDVertex *vertex = esd.GetPrimaryVertexSPD();
1133       Bool_t fVtxOK = kFALSE;
1134       if (vertex->GetNContributors()>0) 
1135       {
1136            TString vtxTyp = vertex->GetTitle();
1137            if ( !vtxTyp.Contains("vertexer: Z") || (vertex->GetDispersion()<0.04 && vertex->GetZRes()<0.25))
1138            {     
1139                 fVtxOK = kTRUE;
1140                 ip.SetX(vertex->GetX());
1141                 ip.SetY(vertex->GetY());
1142                 ip.SetZ(vertex->GetZ());                
1143            }    
1144       }
1145    return fVtxOK;       
1146         
1147 }
1148
1149 //____________________________________________________________________
1150 Bool_t
1151 AliFMDEventInspector::CheckVertex(const AliESDEvent& esd, 
1152                                   TVector3& ip) const
1153 {
1154   // Use standard SPD vertex (perhaps preferable for Pb+Pb)
1155   // Get the vertex 
1156   const AliESDVertex* vertex = esd.GetPrimaryVertexSPD();
1157   if (!vertex) { 
1158     if (fDebug > 2) {
1159       AliWarning("No SPD vertex found in ESD"); }
1160     return false;
1161   }
1162     
1163   // Check that enough tracklets contributed 
1164   if(vertex->GetNContributors() <= 0) {
1165     DMSG(fDebug,2,"Number of contributors to vertex is %d<=0",
1166          vertex->GetNContributors());
1167     ip.SetZ(0);
1168     return false;
1169   } 
1170   // Check that the uncertainty isn't too large 
1171   if (vertex->GetZRes() > fMaxVzErr) { 
1172     DMSG(fDebug,2,"Uncertaintity in Z of vertex is too large %f > %f", 
1173          vertex->GetZRes(), fMaxVzErr);
1174     return false;
1175   }
1176     
1177   // Get the z coordiante 
1178   ip.SetZ(vertex->GetZ());
1179   const AliESDVertex* vertexXY = esd.GetPrimaryVertex();
1180     
1181   
1182   if(!vertexXY->IsFromVertexerZ()) {
1183     ip.SetX(vertexXY->GetX());
1184     ip.SetY(vertexXY->GetY());
1185   }
1186   return true;
1187 }
1188
1189 //____________________________________________________________________
1190 Bool_t
1191 AliFMDEventInspector::ReadRunDetails(const AliESDEvent* esd)
1192 {
1193   // 
1194   // Read the collision system, collision energy, and L3 field setting
1195   // from the ESD
1196   // 
1197   // Parameters:
1198   //   esd ESD to get information from 
1199   // 
1200   // Return:
1201   //    true on success, false 
1202   //
1203   // AliInfo(Form("Parameters from 1st ESD event: cms=%s, sNN=%f, field=%f",
1204   //           esd->GetBeamType(), 2*esd->GetBeamEnergy(), 
1205   //           esd->GetMagneticField()));
1206   DGUARD(fDebug,2,"Read the run details in AliFMDEventInspector");
1207   const char* sys  = esd->GetBeamType();
1208   Float_t     cms  = 2 * esd->GetBeamEnergy();
1209   Float_t     fld  = esd->GetMagneticField();
1210   fCollisionSystem = AliForwardUtil::ParseCollisionSystem(sys);
1211   fEnergy          = AliForwardUtil::ParseCenterOfMassEnergy(fCollisionSystem, 
1212                                                              cms);
1213   fField           = AliForwardUtil::ParseMagneticField(fld);
1214
1215   StoreInformation(esd->GetRunNumber());
1216   if (fCollisionSystem   == AliForwardUtil::kUnknown) { 
1217     AliWarningF("Unknown collision system: %s - please check", sys);
1218     return false;
1219   }
1220   if (fEnergy            <= 0) {
1221     AliWarningF("Unknown CMS energy: %f (%d) - please check", cms, fEnergy);
1222     return false;
1223   }
1224   if (TMath::Abs(fField) >  10) {
1225     AliWarningF("Unknown L3 field setting: %f (%d) - please check", fld,fField);
1226     return false;
1227   }
1228
1229   return true;
1230 }
1231
1232
1233 //____________________________________________________________________
1234 const Char_t*
1235 AliFMDEventInspector::CodeString(UInt_t code)
1236 {
1237   static TString s;
1238   s = "";
1239   if (code & kNoEvent)    s.Append("NOEVENT ");
1240   if (code & kNoTriggers) s.Append("NOTRIGGERS ");
1241   if (code & kNoSPD)      s.Append("NOSPD ");
1242   if (code & kNoFMD)      s.Append("NOFMD ");
1243   if (code & kNoVertex)   s.Append("NOVERTEX ");
1244   if (code & kBadVertex)  s.Append("BADVERTEX ");
1245   return s.Data();
1246 }
1247 //____________________________________________________________________
1248 void
1249 AliFMDEventInspector::Print(Option_t*) const
1250 {
1251   // 
1252   // Print information
1253   // 
1254   //   option Not used 
1255   //
1256   char ind[gROOT->GetDirLevel()+1];
1257   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
1258   ind[gROOT->GetDirLevel()] = '\0';
1259   TString sNN(AliForwardUtil::CenterOfMassEnergyString(fEnergy));
1260   sNN.Strip(TString::kBoth, '0');
1261   sNN.ReplaceAll("GeV", " GeV");
1262   TString field(AliForwardUtil::MagneticFieldString(fField));
1263   field.ReplaceAll("p",  "+");
1264   field.ReplaceAll("m",  "-");
1265   field.ReplaceAll("kG", " kG");
1266   
1267   std::cout << std::boolalpha 
1268             << ind << ClassName() << ": " << GetName() << '\n'
1269             << ind << " Vertex bins:            " << fVtxAxis.GetNbins() << '\n'
1270             << ind << " Vertex range:           [" << fVtxAxis.GetXmin() 
1271             << "," << fVtxAxis.GetXmax() << "]\n"
1272             << ind << " Low flux cut:           " << fLowFluxCut << '\n'
1273             << ind << " Max(delta v_z):         " << fMaxVzErr << " cm\n"
1274             << ind << " Min(nContrib_pileup):   " << fMinPileupContrib << '\n'
1275             << ind << " Min(v-pileup):          " << fMinPileupDistance << '\n'
1276             << ind << " System:                 " 
1277             << AliForwardUtil::CollisionSystemString(fCollisionSystem) << '\n'
1278             << ind << " CMS energy per nucleon: " << sNN << '\n'
1279             << ind << " Field:                  " <<  field << '\n'
1280             << ind << " Satellite events:       " << fUseDisplacedVertices<<'\n'
1281             << ind << " Centrality method:      " << fCentMethod << '\n'
1282             << std::noboolalpha;
1283   if (!fCentAxis) { std::cout << std::flush; return; }
1284   Int_t nBin = fCentAxis->GetNbins();
1285   std::cout << ind << " Centrality axis:        " << nBin << " bins"
1286             << std::flush;
1287   for (Int_t i = 0; i < nBin; i++) { 
1288     if ((i % 10) == 0) std::cout << '\n' << ind << "  ";
1289     std::cout << std::setw(5) << fCentAxis->GetBinLowEdge(i+1) << '-';
1290   }
1291   std::cout << std::setw(5) << fCentAxis->GetBinUpEdge(nBin) << std::endl;
1292 }
1293
1294   
1295 //
1296 // EOF
1297 //
1298