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