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