]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/AliForwardMultiplicityTask.cxx
Merge branch 'workdir'
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / AliForwardMultiplicityTask.cxx
1 // 
2 // Calculate the multiplicity in the forward regions event-by-event 
3 // 
4 // Inputs: 
5 //   - AliESDEvent 
6 //
7 // Outputs: 
8 //   - AliAODForwardMult 
9 // 
10 // Histograms 
11 //   
12 // Corrections used 
13 //
14 #include "AliForwardMultiplicityTask.h"
15 #include "AliTriggerAnalysis.h"
16 #include "AliPhysicsSelection.h"
17 #include "AliLog.h"
18 #include "AliESDEvent.h"
19 #include "AliAODHandler.h"
20 #include "AliMultiplicity.h"
21 #include "AliInputEventHandler.h"
22 #include "AliForwardCorrectionManager.h"
23 #include "AliAnalysisManager.h"
24 #include <TH1.h>
25 #include <TH3D.h>
26 #include <TDirectory.h>
27 #include <TTree.h>
28 #include <TROOT.h>
29 #include <TStopwatch.h>
30 #include <TProfile.h>
31
32 //====================================================================
33 AliForwardMultiplicityTask::AliForwardMultiplicityTask()
34   : AliForwardMultiplicityBase(),
35     fESDFMD(),
36     fEventInspector(),
37     fSharingFilter(),
38     fDensityCalculator(),
39     fCorrections(),
40     fHistCollector(),
41     fEventPlaneFinder()
42 {
43   // 
44   // Constructor
45   //
46   DGUARD(fDebug, 3,"Default CTOR of AliForwardMultiplicityTask");
47 }
48
49 //____________________________________________________________________
50 AliForwardMultiplicityTask::AliForwardMultiplicityTask(const char* name)
51   : AliForwardMultiplicityBase(name),
52     fESDFMD(),
53     fEventInspector("event"),
54     fSharingFilter("sharing"), 
55     fDensityCalculator("density"),
56     fCorrections("corrections"),
57     fHistCollector("collector"),
58     fEventPlaneFinder("eventplane")
59 {
60   // 
61   // Constructor 
62   // 
63   // Parameters:
64   //    name Name of task 
65   //
66   DGUARD(fDebug, 3,"named CTOR of AliForwardMultiplicityTask: %s", name);
67 }
68
69
70 //____________________________________________________________________
71 Bool_t
72 AliForwardMultiplicityTask::PreEvent()
73 {
74   // Clear stuff 
75   fHistos.Clear();
76   fESDFMD.Clear();
77   fAODFMD.Clear();
78   fAODEP.Clear();
79   return true;
80 }
81 //____________________________________________________________________
82 Bool_t
83 AliForwardMultiplicityTask::Event(AliESDEvent& esd)
84 {
85   // 
86   // Process each event 
87   // 
88   // Parameters:
89   //    option Not used
90   //  
91   TStopwatch total;
92   TStopwatch individual;
93   if (fDoTiming) total.Start(true);
94   
95   DGUARD(fDebug,1,"Process the input event");
96
97   // Inspect the event
98   if (fDoTiming) individual.Start(true);
99   Bool_t   lowFlux   = kFALSE;
100   UInt_t   triggers  = 0;
101   UShort_t ivz       = 0;
102   TVector3 ip;
103   Double_t cent      = -1;
104   UShort_t nClusters = 0;
105   UInt_t   found     = fEventInspector.Process(&esd, triggers, lowFlux, 
106                                                ivz, ip, cent, nClusters);
107   if (fDoTiming) fHTiming->Fill(kTimingEventInspector, individual.CpuTime());
108   
109   if (found & AliFMDEventInspector::kNoEvent)    return false;
110   if (found & AliFMDEventInspector::kNoTriggers) return false;
111
112   // Set trigger bits, and mark this event for storage 
113   fAODFMD.SetTriggerBits(triggers);
114   fAODFMD.SetSNN(fEventInspector.GetEnergy());
115   fAODFMD.SetSystem(fEventInspector.GetCollisionSystem());
116   fAODFMD.SetCentrality(cent);
117   fAODFMD.SetNClusters(nClusters);
118   MarkEventForStore();
119  
120   // Do not check if SPD data is there - potential bias 
121   // if (found & AliFMDEventInspector::kNoSPD)      return false;
122   if (found    & AliFMDEventInspector::kNoFMD)      return false;
123   if (found    & AliFMDEventInspector::kNoVertex)   return false;
124   if (triggers & AliAODForwardMult::kPileUp)        return false;
125   fAODFMD.SetIpZ(ip.Z());
126   if (found & AliFMDEventInspector::kBadVertex)     return false;
127
128   // We we do not want to use low flux specific code, we disable it here. 
129   if (!fEnableLowFlux) lowFlux = false;
130
131   // Get FMD data 
132   AliESDFMD* esdFMD = esd.GetFMDData();  
133
134   // Apply the sharing filter (or hit merging or clustering if you like)
135   if (fDoTiming) individual.Start(true);
136   if (!fSharingFilter.Filter(*esdFMD, lowFlux, fESDFMD, ip.Z())) { 
137     AliWarning("Sharing filter failed!");
138     return false;
139   }
140   if (fDoTiming) fHTiming->Fill(kTimingSharingFilter, individual.CpuTime());
141   
142   // Calculate the inclusive charged particle density 
143   if (fDoTiming) individual.Start(true);
144   if (!fDensityCalculator.Calculate(fESDFMD, fHistos, lowFlux, cent, ip)) { 
145     // if (!fDensityCalculator.Calculate(*esdFMD, fHistos, ivz, lowFlux)) { 
146     AliWarning("Density calculator failed!");
147     return false;
148   }
149   if (fDoTiming) fHTiming->Fill(kTimingDensityCalculator,individual.CpuTime());
150
151   // Check if we should do the event plane finder
152   if (fEventInspector.GetCollisionSystem() == AliFMDEventInspector::kPbPb) {
153     if (fDoTiming) individual.Start(true);
154     if (!fEventPlaneFinder.FindEventplane(&esd, fAODEP, 
155                                           &(fAODFMD.GetHistogram()), &fHistos))
156       AliWarning("Eventplane finder failed!");
157     if (fDoTiming) fHTiming->Fill(kTimingEventPlaneFinder,individual.CpuTime());
158   }
159   
160   // Check how many rings have been marked for skipping 
161   Int_t nSkip = 0;
162   for (UShort_t d=1; d<=3; d++) { 
163     for (UShort_t q=0; q<=(d/2); q++) { 
164       TH2D* h = fHistos.Get(d,q == 0 ? 'I' : 'O');
165       if (h && h->TestBit(AliForwardUtil::kSkipRing)) nSkip++;
166     }
167   }
168   if (nSkip > 0) 
169     // Skip the rest if we have too many outliers 
170     return false;
171   
172   // Do the secondary and other corrections. 
173   if (fDoTiming) individual.Start(true);
174   if (!fCorrections.Correct(fHistos, ivz)) { 
175     AliWarning("Corrections failed");
176     return false;
177   }
178   if (fDoTiming) fHTiming->Fill(kTimingCorrections, individual.CpuTime());
179
180   // Collect our `super' histogram 
181   if (fDoTiming) individual.Start(true);
182   if (!fHistCollector.Collect(fHistos, fRingSums, 
183                               ivz, fAODFMD.GetHistogram(),
184                               fAODFMD.GetCentrality())) {
185     AliWarning("Histogram collector failed");
186     return false;
187   }
188   if (fDoTiming) fHTiming->Fill(kTimingHistCollector, individual.CpuTime());
189
190   if (fAODFMD.IsTriggerBits(AliAODForwardMult::kInel) && nSkip < 1) 
191     fHData->Add(&(fAODFMD.GetHistogram()));
192
193   if (fDoTiming) fHTiming->Fill(kTimingTotal, total.CpuTime());
194   
195   return true;
196 }
197
198
199 //
200 // EOF
201 //