]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/AliFMDEventInspector.cxx
AliFMDEventInspector, AliFMDDensityCalculator, AliFMDSharingFilter:
[u/mrichter/AliRoot.git] / PWG2 / 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 "AliInputEventHandler.h"
23 #include "AliTriggerAnalysis.h"
24 #include "AliPhysicsSelection.h"
25 #include "AliAODForwardMult.h"
26 #include "AliForwardUtil.h"
27 #include "AliCentrality.h"
28 #include <TH1.h>
29 #include <TList.h>
30 #include <TDirectory.h>
31 #include <TROOT.h>
32 #include <iostream>
33 #include <iomanip>
34
35 //====================================================================
36 AliFMDEventInspector::AliFMDEventInspector()
37   : TNamed(),
38     fHEventsTr(0), 
39     fHEventsTrVtx(0),
40     fHEventsAccepted(0),
41     fHTriggers(0),
42     fHType(0),
43     fHWords(0),
44     fHCent(0),
45     fHCentVsQual(0),
46     fLowFluxCut(1000),
47     fMaxVzErr(0.2),
48     fList(0),
49     fEnergy(0),
50     fField(999), 
51     fCollisionSystem(kUnknown),
52     fDebug(0),
53     fCentAxis(0),
54     fVtxAxis(10,-10,10)
55 {
56   // 
57   // Constructor 
58   //
59 }
60
61 //____________________________________________________________________
62 AliFMDEventInspector::AliFMDEventInspector(const char* name)
63   : TNamed("fmdEventInspector", name),
64     fHEventsTr(0), 
65     fHEventsTrVtx(0), 
66     fHEventsAccepted(0),
67     fHTriggers(0),
68     fHType(0),
69     fHWords(0),
70     fHCent(0),
71     fHCentVsQual(0),
72     fLowFluxCut(1000),
73     fMaxVzErr(0.2),
74     fList(0),
75     fEnergy(0),
76     fField(999), 
77     fCollisionSystem(kUnknown),
78     fDebug(0),
79     fCentAxis(0),
80     fVtxAxis(10,-10,10)
81 {
82   // 
83   // Constructor 
84   // 
85   // Parameters:
86   //   name Name of object
87   //
88 }
89
90 //____________________________________________________________________
91 AliFMDEventInspector::AliFMDEventInspector(const AliFMDEventInspector& o)
92   : TNamed(o), 
93     fHEventsTr(o.fHEventsTr), 
94     fHEventsTrVtx(o.fHEventsTrVtx), 
95     fHEventsAccepted(o.fHEventsAccepted),
96     fHTriggers(o.fHTriggers),
97     fHType(o.fHType),
98     fHWords(o.fHWords),
99     fHCent(o.fHCent),
100     fHCentVsQual(o.fHCentVsQual),
101     fLowFluxCut(o.fLowFluxCut),
102     fMaxVzErr(o.fMaxVzErr),
103     fList(o.fList),
104     fEnergy(o.fEnergy),
105     fField(o.fField), 
106     fCollisionSystem(o.fCollisionSystem),
107     fDebug(0),
108     fCentAxis(0),
109     fVtxAxis(o.fVtxAxis)
110 {
111   // 
112   // Copy constructor 
113   // 
114   // Parameters:
115   //   o Object to copy from 
116   //
117 }
118
119 //____________________________________________________________________
120 AliFMDEventInspector::~AliFMDEventInspector()
121 {
122   // 
123   // Destructor 
124   //
125   if (fList)         delete fList;
126 }
127 //____________________________________________________________________
128 AliFMDEventInspector&
129 AliFMDEventInspector::operator=(const AliFMDEventInspector& o)
130 {
131   // 
132   // Assignement operator
133   // 
134   // Parameters:
135   //   o Object to assign from 
136   // 
137   // Return:
138   //    Reference to this object
139   //
140   TNamed::operator=(o);
141   fHEventsTr         = o.fHEventsTr;
142   fHEventsTrVtx      = o.fHEventsTrVtx;
143   fHEventsAccepted   = o.fHEventsAccepted;
144   fHTriggers         = o.fHTriggers;
145   fHType             = o.fHType;
146   fHWords            = o.fHWords;
147   fHCent             = o.fHCent;
148   fHCentVsQual       = o.fHCentVsQual;
149   fLowFluxCut        = o.fLowFluxCut;
150   fMaxVzErr          = o.fMaxVzErr;
151   fDebug             = o.fDebug;
152   fList              = (o.fList ? new TList : 0);
153   fEnergy            = o.fEnergy;
154   fField             = o.fField;
155   fCollisionSystem   = o.fCollisionSystem;
156   fVtxAxis.Set(o.fVtxAxis.GetNbins(), o.fVtxAxis.GetXmin(), 
157                o.fVtxAxis.GetXmax());
158   if (fList) { 
159     fList->SetName(GetName());
160     if (fHEventsTr)    fList->Add(fHEventsTr);
161     if (fHEventsTrVtx) fList->Add(fHEventsTrVtx);
162     if (fHTriggers)    fList->Add(fHTriggers);
163     if (fHType)        fList->Add(fHType);
164     if (fHWords)       fList->Add(fHWords);
165     if (fHCent)        fList->Add(fHCent);
166     if (fHCentVsQual)  fList->Add(fHCentVsQual);
167   }
168   return *this;
169 }
170
171 //____________________________________________________________________
172 Bool_t
173 AliFMDEventInspector::FetchHistograms(const TList* d, 
174                                       TH1I*& hEventsTr, 
175                                       TH1I*& hEventsTrVtx, 
176                                       TH1I*& hTriggers) const
177 {
178   // 
179   // Fetch our histograms from the passed list 
180   // 
181   // Parameters:
182   //   d             Input
183   //   hEventsTr     On return, pointer to histogram, or null
184   //   hEventsTrVtx  On return, pointer to histogram, or null
185   //   hTriggers     On return, pointer to histogram, or null
186   // 
187   // Return:
188   //    true on success, false otherwise 
189   //
190   hEventsTr    = 0;
191   hEventsTrVtx = 0;
192   hTriggers    = 0;
193   TList* dd    = dynamic_cast<TList*>(d->FindObject(GetName()));
194   if (!dd) return kFALSE;
195   
196   hEventsTr    = dynamic_cast<TH1I*>(dd->FindObject("nEventsTr"));
197   hEventsTrVtx = dynamic_cast<TH1I*>(dd->FindObject("nEventsTrVtx"));
198   hTriggers    = dynamic_cast<TH1I*>(dd->FindObject("triggers"));
199
200   if (!hEventsTr || !hEventsTrVtx || !hTriggers) return kFALSE;
201   return kTRUE;
202 }
203 //____________________________________________________________________
204 void
205 AliFMDEventInspector::Init(const TAxis& vtxAxis)
206 {
207   // 
208   // Initialize the object 
209   // 
210   // Parameters:
211   //   vtxAxis Vertex axis in use 
212   //
213   
214   // -1.5 -0.5 0.5 1.5 ... 89.5 ... 100.5
215   // ----- 92 number --------- ---- 1 ---
216   TArrayD limits(93);
217   for (Int_t i = 0; i < 92; i++) limits[i] = -1.5 + i;
218   limits[92] = 100.5;
219
220   fVtxAxis.Set(vtxAxis.GetNbins(), vtxAxis.GetXmin(), vtxAxis.GetXmax());
221   
222   fCentAxis  = new TAxis(limits.GetSize()-1, limits.GetArray());
223   fHEventsTr = new TH1I("nEventsTr", "Number of events w/trigger", 
224                         4*vtxAxis.GetNbins(), 
225                         2*vtxAxis.GetXmin(), 
226                         2*vtxAxis.GetXmax());
227   fHEventsTr->SetXTitle("v_{z} [cm]");
228   fHEventsTr->SetYTitle("# of events");
229   fHEventsTr->SetFillColor(kRed+1);
230   fHEventsTr->SetFillStyle(3001);
231   fHEventsTr->SetDirectory(0);
232   // fHEventsTr->Sumw2();
233   fList->Add(fHEventsTr);
234
235   fHEventsTrVtx = static_cast<TH1I*>(fHEventsTr->Clone("nEventsTrVtx"));
236   fHEventsTrVtx->SetTitle("Number of events w/trigger and vertex"); 
237   fHEventsTrVtx->SetFillColor(kBlue+1);
238   fHEventsTrVtx->SetDirectory(0);
239   // fHEventsTrVtx->Sumw2();
240   fList->Add(fHEventsTrVtx);
241
242   fHEventsAccepted = new TH1I("nEventsAccepted", 
243                               "Number of events  w/trigger and vertex in range",
244                               2*vtxAxis.GetNbins(), 
245                               2*vtxAxis.GetXmin(), 
246                               2*vtxAxis.GetXmax());
247   fHEventsAccepted->SetXTitle("v_{z} [cm]");
248   fHEventsAccepted->SetYTitle("# of events");
249   fHEventsAccepted->SetFillColor(kGreen+1);
250   fHEventsAccepted->SetFillStyle(3001);
251   fHEventsAccepted->SetDirectory(0);
252   // fHEventsAccepted->Sumw2();
253   fList->Add(fHEventsAccepted);
254                               
255       
256   fHTriggers = new TH1I("triggers", "Triggers", kOffline+1, 0, kOffline+1);
257   fHTriggers->SetFillColor(kRed+1);
258   fHTriggers->SetFillStyle(3001);
259   fHTriggers->SetStats(0);
260   fHTriggers->SetDirectory(0);
261   fHTriggers->GetXaxis()->SetBinLabel(kInel   +1,"INEL");
262   fHTriggers->GetXaxis()->SetBinLabel(kInelGt0+1,"INEL>0");
263   fHTriggers->GetXaxis()->SetBinLabel(kNSD    +1,"NSD");
264   fHTriggers->GetXaxis()->SetBinLabel(kEmpty  +1,"Empty");
265   fHTriggers->GetXaxis()->SetBinLabel(kA      +1,"A");
266   fHTriggers->GetXaxis()->SetBinLabel(kB      +1,"B");
267   fHTriggers->GetXaxis()->SetBinLabel(kC      +1,"C");
268   fHTriggers->GetXaxis()->SetBinLabel(kE      +1,"E");
269   fHTriggers->GetXaxis()->SetBinLabel(kPileUp +1,"Pileup");
270   fHTriggers->GetXaxis()->SetBinLabel(kMCNSD  +1,"NSD_{MC}");
271   fHTriggers->GetXaxis()->SetBinLabel(kOffline+1,"Offline");
272   fList->Add(fHTriggers);
273
274   fHType = new TH1I("type", Form("Event type (cut: SPD mult>%d)", 
275                                  fLowFluxCut), 2, -.5, 1.5);
276   fHType->SetFillColor(kRed+1);
277   fHType->SetFillStyle(3001);
278   fHType->SetStats(0);
279   fHType->SetDirectory(0);
280   fHType->GetXaxis()->SetBinLabel(1,"Low-flux");
281   fHType->GetXaxis()->SetBinLabel(2,"High-flux");
282   fList->Add(fHType);
283
284
285   fHWords = new TH1I("words", "Trigger words seen", 1, 0, 0); 
286   fHWords->SetFillColor(kBlue+1);
287   fHWords->SetFillStyle(3001);
288   fHWords->SetStats(0);
289   fHWords->SetDirectory(0);
290   fHWords->SetBit(TH1::kCanRebin);
291   fList->Add(fHWords);
292
293   fHCent = new TH1F("cent", "Centrality", limits.GetSize()-1,limits.GetArray());
294   fHCent->SetFillColor(kBlue+1);
295   fHCent->SetFillStyle(3001);
296   fHCent->SetStats(0);
297   fHCent->SetDirectory(0);
298   fHCent->SetXTitle("Centrality [%]");
299   fHCent->SetYTitle("Events");
300   fList->Add(fHCent);
301
302   fHCentVsQual = new TH2F("centVsQuality", "Quality vs Centrality", 
303                           5, 0, 5, limits.GetSize()-1, limits.GetArray());
304   fHCentVsQual->SetXTitle("Quality");
305   fHCentVsQual->SetYTitle("Centrality [%]");
306   fHCentVsQual->SetZTitle("Events");
307   fHCentVsQual->GetXaxis()->SetBinLabel(1, "OK");
308   fHCentVsQual->GetXaxis()->SetBinLabel(2, "Outside v_{z} cut");
309   fHCentVsQual->GetXaxis()->SetBinLabel(3, "V0 vs SPD outlier");
310   fHCentVsQual->GetXaxis()->SetBinLabel(4, "V0 vs TPC outlier");
311   fHCentVsQual->GetXaxis()->SetBinLabel(5, "V0 vs ZDC outlier");
312   fList->Add(fHCentVsQual);
313 }
314
315 //____________________________________________________________________
316 void
317 AliFMDEventInspector::StoreInformation()
318 {
319   // Write TNamed objects to output list containing information about
320   // the running conditions 
321   if (!fList) return;
322
323   TNamed* sys = new TNamed("sys", "");
324   TNamed* sNN = new TNamed("sNN", "");
325   TNamed* fld = new TNamed("field", "");
326   sys->SetTitle(AliForwardUtil::CollisionSystemString(fCollisionSystem));
327   sNN->SetTitle(AliForwardUtil::CenterOfMassEnergyString(fEnergy));
328   fld->SetTitle(AliForwardUtil::MagneticFieldString(fField));
329   sys->SetUniqueID(fCollisionSystem);
330   sNN->SetUniqueID(fEnergy);
331   fld->SetUniqueID(fField);
332
333   fList->Add(sys);
334   fList->Add(sNN);
335   fList->Add(fld);
336 }
337
338 //____________________________________________________________________
339 void
340 AliFMDEventInspector::DefineOutput(TList* dir)
341 {
342   // 
343   // Define the output histograms.  These are put in a sub list of the
344   // passed list.   The histograms are merged before the parent task calls 
345   // AliAnalysisTaskSE::Terminate 
346   // 
347   //   dir Directory to add to 
348   //
349   fList = new TList;
350   fList->SetName(GetName());
351   dir->Add(fList);
352 }
353
354 //____________________________________________________________________
355 UInt_t
356 AliFMDEventInspector::Process(const AliESDEvent* event, 
357                               UInt_t&            triggers,
358                               Bool_t&            lowFlux,
359                               UShort_t&          ivz, 
360                               Double_t&          vz,
361                               Double_t&          cent,
362                               UShort_t&          nClusters)
363 {
364   // 
365   // Process the event 
366   // 
367   // Parameters:
368   //   event     Input event 
369   //   triggers  On return, the triggers fired 
370   //   lowFlux   On return, true if the event is considered a low-flux 
371   //                  event (according to the setting of fLowFluxCut) 
372   //   ivz       On return, the found vertex bin (1-based).  A zero
373   //                  means outside of the defined vertex range
374   //   vz        On return, the z position of the interaction
375   //   cent      On return, the centrality - if not available < 0
376   // 
377   // Return:
378   //    0 (or kOk) on success, otherwise a bit mask of error codes 
379   //
380
381   // --- Check that we have an event ---------------------------------
382   if (!event) { 
383     AliWarning("No ESD event found for input event");
384     return kNoEvent;
385   }
386
387   // --- Read trigger information from the ESD and store in AOD object
388   if (!ReadTriggers(event, triggers, nClusters)) { 
389     if (fDebug > 2) {
390       AliWarning("Failed to read triggers from ESD"); }
391     return kNoTriggers;
392   }
393
394   // --- Check if this is a high-flux event --------------------------
395   const AliMultiplicity* testmult = event->GetMultiplicity();
396   if (!testmult) {
397     if (fDebug > 3) {
398       AliWarning("No central multiplicity object found"); }
399   }
400   else 
401     lowFlux = testmult->GetNumberOfTracklets() < fLowFluxCut;
402
403   fHType->Fill(lowFlux ? 0 : 1);
404   
405   // --- Read centrality information 
406   cent          = -10;
407   UShort_t qual = 0;
408   if (!ReadCentrality(event, cent, qual)) {
409     if (fDebug > 3) 
410       AliWarning("Failed to get centrality");
411   }
412   fHCent->Fill(cent);
413   if (qual == 0) fHCentVsQual->Fill(0., cent);
414   else { 
415     for (UShort_t i = 0; i < 4; i++) 
416       if (qual & (1 << i)) fHCentVsQual->Fill(Double_t(i+1), cent);
417   }
418
419   // --- Get the vertex information ----------------------------------
420   vz          = 0;
421   Bool_t vzOk = ReadVertex(event, vz);
422
423   fHEventsTr->Fill(vz);
424   if (!vzOk) { 
425     if (fDebug > 3) {
426       AliWarning("Failed to read vertex from ESD"); }
427     return kNoVertex;
428   }
429   fHEventsTrVtx->Fill(vz);
430
431   // --- Get the vertex bin ------------------------------------------
432   ivz = fVtxAxis.FindBin(vz);
433   if (ivz <= 0 || ivz > fVtxAxis.GetNbins()) { 
434     if (fDebug > 3) {
435       AliWarning(Form("Vertex @ %f outside of range [%f,%f]", 
436                       vz, fVtxAxis.GetXmin(), fVtxAxis.GetXmax())); 
437     }
438     ivz = 0;
439     return kBadVertex;
440   }
441   fHEventsAccepted->Fill(vz);
442   
443   // --- Check the FMD ESD data --------------------------------------
444   if (!event->GetFMDData()) { 
445     if (fDebug > 3) {
446       AliWarning("No FMD data found in ESD"); }
447     return kNoFMD;
448   }
449
450   
451   return kOk;
452 }
453
454 //____________________________________________________________________
455 Bool_t
456 AliFMDEventInspector::ReadCentrality(const AliESDEvent* esd, 
457                                      Double_t& cent, 
458                                      UShort_t& qual) const
459 {
460   // 
461   // Read centrality from event 
462   // 
463   // Parameters:
464   //    esd  Event 
465   //    cent On return, the centrality or negative if not found
466   // 
467   // Return:
468   //    False on error, true otherwise 
469   //
470   cent = -1;
471   qual = 0;
472   AliCentrality* centObj = const_cast<AliESDEvent*>(esd)->GetCentrality();
473   if (!centObj)  return true;
474
475   // AliInfo(Form("Got centrality object %p with quality %d", 
476   //              centObj, centObj->GetQuality()));
477   // centObj->Print();
478   cent = centObj->GetCentralityPercentile("V0M");  
479   qual = centObj->GetQuality();
480
481   return true;
482 }
483
484 //____________________________________________________________________
485 Bool_t
486 AliFMDEventInspector::ReadTriggers(const AliESDEvent* esd, UInt_t& triggers,
487                                    UShort_t& nClusters)
488 {
489   // 
490   // Read the trigger information from the ESD event 
491   // 
492   // Parameters:
493   //   esd        ESD event 
494   //   triggers   On return, contains the trigger bits 
495   // 
496   // Return:
497   //    @c true on success, @c false otherwise 
498   //
499   triggers = 0;
500
501   // Get the analysis manager - should always be there 
502   AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
503   if (!am) { 
504     AliWarning("No analysis manager defined!");
505     return kFALSE;
506   }
507
508   // Get the input handler - should always be there 
509   AliInputEventHandler* ih = 
510     static_cast<AliInputEventHandler*>(am->GetInputEventHandler());
511   if (!ih) { 
512     AliWarning("No input handler");
513     return kFALSE;
514   }
515     
516   // Check if this is a collision candidate (MB)
517   // Note, that we should use the value cached in the input 
518   // handler rather than calling IsCollisionCandiate directly 
519   // on the AliPhysicsSelection obejct.  If we called the latter
520   // then the AliPhysicsSelection object would overcount by a 
521   // factor of 2! :-(
522   Bool_t offline = ih->IsEventSelected();
523   nClusters = 0;
524   if (offline) {
525     triggers |= AliAODForwardMult::kOffline;
526     triggers |= AliAODForwardMult::kInel;
527     fHTriggers->Fill(kOffline+0.5);
528
529     // If this is inel, see if we have a tracklet 
530     const AliMultiplicity* spdmult = esd->GetMultiplicity();
531     if (!spdmult) {
532       AliWarning("No SPD multiplicity");
533     }
534     else { 
535       // Check if we have one or more tracklets 
536       // in the range -1 < eta < 1 to set the INEL>0 
537       // trigger flag. 
538       // 
539       // Also count tracklets as a single cluster 
540       Int_t n = spdmult->GetNumberOfTracklets();
541       for (Int_t j = 0; j < n; j++) { 
542         if(TMath::Abs(spdmult->GetEta(j)) < 1) { 
543           triggers |= AliAODForwardMult::kInelGt0;
544           nClusters++;
545         }
546       }
547       n = spdmult->GetNumberOfSingleClusters();
548       for (Int_t j = 0; j < n; j++) { 
549         Double_t eta = -TMath::Log(TMath::Tan(spdmult->GetThetaSingle(j)/2.));
550         if (TMath::Abs(eta) < 1) nClusters++;
551       }
552     }
553     if (nClusters > 0) triggers |= AliAODForwardMult::kNClusterGt0;
554   }
555
556   // Analyse some trigger stuff 
557   AliTriggerAnalysis ta;
558   if (ta.IsOfflineTriggerFired(esd, AliTriggerAnalysis::kNSD1)) 
559     triggers |= AliAODForwardMult::kNSD;
560  
561     
562   // Check for multiple vertices (pile-up) with at least 3
563   // contributors and at least 0.8cm from the primary vertex
564   Bool_t pileup =  esd->IsPileupFromSPD(3,0.8);
565   if (pileup) {
566     triggers |= AliAODForwardMult::kPileUp;
567     fHTriggers->Fill(kPileUp+.5);
568   }
569
570     
571   // Get trigger stuff 
572   TString trigStr = esd->GetFiredTriggerClasses();
573   // AliWarning(Form("Fired trigger classes: %s", trigStr.Data()));
574   fHWords->Fill(trigStr.Data(), 1);
575 #if 0
576   if (trigStr.Contains("MB1") || trigStr.Contains("MBBG3"))
577     triggers |= AliAOODForwardMult::kB;
578   if (trigStr.Contains("COTA")) 
579     triggers |= AliAODForwardMult::kA;
580   if (trigStr.Contains("COTC")) 
581     triggers |= AliAODForwardMult::kC;
582 #endif
583   if (trigStr.Contains("CBEAMB-ABCE-NOPF-ALL")) {
584     triggers |= AliAODForwardMult::kEmpty;
585     fHTriggers->Fill(kEmpty+.5);
586   }
587
588   // Check for B triggers
589   if (trigStr.Contains("CINT1B-ABCE-NOPF-ALL")   ||   // Early pp
590       trigStr.Contains("CINT1-B-NOPF-ALLNOTRD")  ||   // Late pp 
591       trigStr.Contains("CINT1-B-NOPF-FASTNOTRD") ||   // Late pp 
592       trigStr.Contains("CSMBB-ABCE-NOPF-ALL")    ||   // pp
593       trigStr.Contains("CMBACS2-B-NOPF-ALL")     ||   // PbPb
594       // trigStr.Contains("C0SMH-B-NOPF-ALL")    ||   // PbPb - high mult
595       trigStr.Contains("CMBS2A-B-NOPF-ALL")      ||   // PbPb
596       trigStr.Contains("CMBS2C-B-NOPF-ALL")      ||   // PbPb
597       trigStr.Contains("CMBAC-B-NOPF-ALL")       ||   // PbPb
598       // trigStr.Contains("C0SMH-B-NOPF-ALL")    ||   // PbPb - high mult
599       trigStr.Contains("CMBACS2-B-NOPF-ALLNOTRD")     // PbPb
600       // trigStr.Contains("C0SMH-B-NOPF-ALLNOTRD")    // PbPb - high mult
601       ) {
602     triggers |= AliAODForwardMult::kB;
603     fHTriggers->Fill(kB+.5);
604   }
605   
606   // Check for A triggers
607   if (trigStr.Contains("CINT1A-ABCE-NOPF-ALL")   ||   // Early pp
608       trigStr.Contains("CINT1-AC_NOPF-ALLNOTRD") ||   // Late pp
609       trigStr.Contains("CINT1-AC_NOPF-FASTNOTRD")||   // Late pp
610       (trigStr.Contains("CSMBA-ABCE-NOPF-ALL") && 
611        !(triggers & AliAODForwardMult::kB))      ||   // pp
612       trigStr.Contains("CMBACS2-A-NOPF-ALL")     ||   // PbPb
613       // trigStr.Contains("C0SMH-A-NOPF-ALL")    ||   // PbPb - high mult
614       trigStr.Contains("CMBS2A-A-NOPF-ALL")      ||   // PbPb
615       trigStr.Contains("CMBS2C-A-NOPF-ALL")      ||   // PbPb
616       trigStr.Contains("CMBAC-A-NOPF-ALL")       ||   // PbPb
617       // trigStr.Contains("C0SMH-A-NOPF-ALL")    ||   // PbPb - high mult
618       trigStr.Contains("CMBACS2-A-NOPF-ALLNOTRD")     // PbPb
619       // trigStr.Contains("C0SMH-A-NOPF-ALLNOTRD")    // PbPb - high mult
620       ) {
621     triggers |= AliAODForwardMult::kA;
622     fHTriggers->Fill(kA+.5);
623   }
624
625   // Check for C triggers
626   if (trigStr.Contains("CINT1C-ABCE-NOPF-ALL")   ||  // Early pp
627       (trigStr.Contains("CSMBC-ABCE-NOPF-ALL") && 
628        !(triggers & AliAODForwardMult::kB))      ||   // pp
629       trigStr.Contains("CMBACS2-C-NOPF-ALL")     ||   // PbPb
630       // trigStr.Contains("C0SMH-B-NOPF-ALL")    ||   // PbPb - high mult
631       trigStr.Contains("CMBS2A-C-NOPF-ALL")      ||   // PbPb
632       trigStr.Contains("CMBS2C-C-NOPF-ALL")      ||   // PbPb
633       trigStr.Contains("CMBAC-C-NOPF-ALL")       ||   // PbPb
634       // trigStr.Contains("C0SMH-B-NOPF-ALL")    ||   // PbPb - high mult
635       trigStr.Contains("CMBACS2-C-NOPF-ALLNOTRD")     // PbPb
636       // trigStr.Contains("C0SMH-B-NOPF-ALLNOTRD")    // PbPb - high mult
637       ) {
638     triggers |= AliAODForwardMult::kC;
639     fHTriggers->Fill(kC+.5);
640   }
641
642   // Check for E triggers 
643   if (trigStr.Contains("CINT1-E-NOPF-ALL")       ||   // Early pp 
644       trigStr.Contains("CINT1-E-NOPF-ALLNOTRD")  ||   // Late pp 
645       trigStr.Contains("CINT1-E-NOPF-FASTNOTRD") ||   // Late pp 
646       trigStr.Contains("CMBACS2-E-NOPF-ALL")     ||   // PbPb
647       // trigStr.Contains("C0SMH-B-NOPF-ALL")    ||   // PbPb - high mult
648       trigStr.Contains("CMBS2A-E-NOPF-ALL")      ||   // PbPb
649       trigStr.Contains("CMBS2C-E-NOPF-ALL")      ||   // PbPb
650       trigStr.Contains("CMBAC-E-NOPF-ALL")       ||   // PbPb
651       // trigStr.Contains("C0SMH-B-NOPF-ALL")    ||   // PbPb - high mult
652       trigStr.Contains("CMBACS2-E-NOPF-ALLNOTRD")     // PbPb
653       // trigStr.Contains("C0SMH-B-NOPF-ALLNOTRD")    // PbPb - high mult
654       ) {
655     triggers |= AliAODForwardMult::kE;
656     fHTriggers->Fill(kE+.5);
657   }
658
659   // Now check - if we have a collision - for offline triggers and
660   // fill histogram.
661   if (triggers & AliAODForwardMult::kB) {
662     if (triggers & AliAODForwardMult::kInel) 
663       fHTriggers->Fill(kInel);
664     
665     if (triggers & AliAODForwardMult::kInelGt0)
666       fHTriggers->Fill(kInelGt0+.5);
667     
668     if (triggers & AliAODForwardMult::kNSD)
669       fHTriggers->Fill(kNSD+.5);
670   }
671   
672   return kTRUE;
673 }
674 //____________________________________________________________________
675 Bool_t
676 AliFMDEventInspector::ReadVertex(const AliESDEvent* esd, Double_t& vz)
677 {
678   // 
679   // Read the vertex information from the ESD event 
680   // 
681   // Parameters:
682   //   esd  ESD event 
683   //   vz   On return, the vertex Z position 
684   // 
685   // Return:
686   //    @c true on success, @c false otherwise 
687   //
688   vz = 0;
689 #if 1
690   // This is the code used by the 1st physics people 
691   const AliESDVertex* vertex    = esd->GetPrimaryVertex();
692   if (!vertex  || !vertex->GetStatus()) {
693     if (fDebug > 2) {
694       AliWarning(Form("No primary vertex (%p) or bad status %d", 
695                       vertex, (vertex ? vertex->GetStatus() : -1)));
696     }
697     return false;
698   }
699   const AliESDVertex* vertexSPD = esd->GetPrimaryVertexSPD();
700   if (!vertexSPD || !vertexSPD->GetStatus()) {
701     if (fDebug > 2) {
702       AliWarning(Form("No primary SPD vertex (%p) or bad status %d", 
703                       vertexSPD, (vertexSPD ? vertexSPD->GetStatus() : -1)));
704     }
705     return false;
706   }
707
708   // if vertex is from SPD vertexZ, require more stringent cuts 
709   if (vertex->IsFromVertexerZ()) { 
710     if (vertex->GetDispersion() > fMaxVzErr || 
711         vertex->GetZRes() > 1.25 * fMaxVzErr) {
712       if (fDebug > 2) {
713         AliWarning(Form("Dispersion %f > %f or resolution %f > %f",
714                         vertex->GetDispersion(), fMaxVzErr,
715                         vertex->GetZRes(), 1.25 * fMaxVzErr)); 
716       }
717       return false;
718     }
719   }
720   vz = vertex->GetZ();
721   return true;
722 #else 
723   // Get the vertex 
724   const AliESDVertex* vertex = esd->GetPrimaryVertexSPD();
725   if (!vertex) { 
726     if (fDebug > 2) {
727       AliWarning("No SPD vertex found in ESD"); }
728     return kFALSE;
729   }
730   
731   // Check that enough tracklets contributed 
732   if(vertex->GetNContributors() <= 0) {
733     if (fDebug > 2) {
734       AliWarning(Form("Number of contributors to vertex is %d<=0",
735                       vertex->GetNContributors())); }
736     vz = 0;
737     return kFALSE;
738   } 
739   // Check that the uncertainty isn't too large 
740   if (vertex->GetZRes() > fMaxVzErr) { 
741     if (fDebug > 2) {
742       AliWarning(Form("Uncertaintity in Z of vertex is too large %f > %f", 
743                       vertex->GetZRes(), fMaxVzErr)); }
744     return kFALSE;
745   }
746   
747   // Get the z coordiante 
748   vz = vertex->GetZ();
749   return kTRUE;
750 #endif 
751 }
752   
753 //____________________________________________________________________
754 Bool_t
755 AliFMDEventInspector::ReadRunDetails(const AliESDEvent* esd)
756 {
757   // 
758   // Read the collision system, collision energy, and L3 field setting
759   // from the ESD
760   // 
761   // Parameters:
762   //   esd ESD to get information from 
763   // 
764   // Return:
765   //    true on success, false 
766   //
767   // AliInfo(Form("Parameters from 1st ESD event: cms=%s, sNN=%f, field=%f",
768   //           esd->GetBeamType(), 2*esd->GetBeamEnergy(), 
769   //           esd->GetMagneticField()));
770   fCollisionSystem = 
771     AliForwardUtil::ParseCollisionSystem(esd->GetBeamType());
772   fEnergy          = 
773     AliForwardUtil::ParseCenterOfMassEnergy(fCollisionSystem,   
774                                             2 * esd->GetBeamEnergy());
775   fField           = 
776     AliForwardUtil::ParseMagneticField(esd->GetMagneticField());
777
778   StoreInformation();
779   if (fCollisionSystem   == AliForwardUtil::kUnknown || 
780       fEnergy            <= 0                        || 
781       TMath::Abs(fField) >  10) 
782     return kFALSE;
783
784   return kTRUE;
785 }
786
787 //____________________________________________________________________
788 void
789 AliFMDEventInspector::Print(Option_t*) const
790 {
791   // 
792   // Print information
793   // 
794   //   option Not used 
795   //
796   char ind[gROOT->GetDirLevel()+1];
797   for (Int_t i = 0; i < gROOT->GetDirLevel(); i++) ind[i] = ' ';
798   ind[gROOT->GetDirLevel()] = '\0';
799   TString sNN(AliForwardUtil::CenterOfMassEnergyString(fEnergy));
800   sNN.Strip(TString::kBoth, '0');
801   sNN.ReplaceAll("GeV", " GeV");
802   TString field(AliForwardUtil::MagneticFieldString(fField));
803   field.ReplaceAll("p",  "+");
804   field.ReplaceAll("m",  "-");
805   field.ReplaceAll("kG", " kG");
806   
807   std::cout << ind << ClassName() << ": " << GetName() << '\n'
808             << ind << " Vertex bins:            " << fVtxAxis.GetNbins() << '\n'
809             << ind << " Vertex range:           [" << fVtxAxis.GetXmin() 
810             << "," << fVtxAxis.GetXmax() << "]\n"
811             << ind << " Low flux cut:           " << fLowFluxCut << '\n'
812             << ind << " Max(delta v_z):         " << fMaxVzErr << " cm\n"
813             << ind << " System:                 " 
814             << AliForwardUtil::CollisionSystemString(fCollisionSystem) << '\n'
815             << ind << " CMS energy per nucleon: " << sNN << '\n'
816             << ind << " Field:                  " <<  field << '\n';
817   if (!fCentAxis) { std::cout << std::flush; return; }
818   Int_t nBin = fCentAxis->GetNbins();
819   std::cout << ind << " Centrality axis:        " << nBin << " bins"
820             << std::flush;
821   for (Int_t i = 0; i < nBin; i++) { 
822     if ((i % 10) == 0) std::cout << '\n' << ind << "  ";
823     std::cout << std::setw(5) << fCentAxis->GetBinLowEdge(i+1) << '-';
824   }
825   std::cout << std::setw(5) << fCentAxis->GetBinUpEdge(nBin) << std::endl;
826 }
827
828   
829 //
830 // EOF
831 //
832