]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/AliForwardMultiplicityBase.cxx
Preliminary work to get centrality in
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / AliForwardMultiplicityBase.cxx
1 //====================================================================
2 // 
3 // Base class for classes that calculate the multiplicity in the
4 // forward regions event-by-event
5 // 
6 // Inputs: 
7 //   - AliESDEvent 
8 //
9 // Outputs: 
10 //   - AliAODForwardMult 
11 // 
12 // Histograms 
13 //   
14 // Corrections used 
15 #include "AliForwardMultiplicityBase.h"
16 #include "AliForwardCorrectionManager.h"
17 #include "AliLog.h"
18 #include "AliAODHandler.h"
19 #include "AliInputEventHandler.h"
20 #include "AliAnalysisManager.h"
21 #include "AliFMDEventInspector.h"
22 #include "AliFMDEnergyFitter.h"
23 #include "AliFMDSharingFilter.h"
24 #include "AliFMDDensityCalculator.h"
25 #include "AliFMDCorrector.h"
26 #include "AliFMDHistCollector.h"
27 #include "AliESDEvent.h"
28 #include <TROOT.h>
29 #include <iostream>
30 #include <iomanip>
31
32 //====================================================================
33 AliForwardMultiplicityBase::AliForwardMultiplicityBase(const char* name) 
34   : AliAnalysisTaskSE(name), 
35     fEnableLowFlux(true), 
36     fFirstEvent(true),
37     fCorrManager(0)
38 {
39   // Set our persistent pointer 
40   fCorrManager = &AliForwardCorrectionManager::Instance();
41 }
42
43 //____________________________________________________________________
44 Bool_t 
45 AliForwardMultiplicityBase::CheckCorrections(UInt_t what) const
46 {
47   // 
48   // Check if all needed corrections are there and accounted for.  If not,
49   // do a Fatal exit 
50   // 
51   // Parameters:
52   //    what Which corrections is needed
53   // 
54   // Return:
55   //    true if all present, false otherwise
56   //  
57
58   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
59   // Check that we have the energy loss fits, needed by 
60   //   AliFMDSharingFilter 
61   //   AliFMDDensityCalculator 
62   if (what & AliForwardCorrectionManager::kELossFits && !fcm.GetELossFit()) { 
63     AliFatal(Form("No energy loss fits"));
64     return false;
65   }
66   // Check that we have the double hit correction - (optionally) used by 
67   //  AliFMDDensityCalculator 
68   if (what & AliForwardCorrectionManager::kDoubleHit && !fcm.GetDoubleHit()) {
69     AliFatal("No double hit corrections"); 
70     return false;
71   }
72   // Check that we have the secondary maps, needed by 
73   //   AliFMDCorrector 
74   //   AliFMDHistCollector
75   if (what & AliForwardCorrectionManager::kSecondaryMap && 
76       !fcm.GetSecondaryMap()) {
77     AliFatal("No secondary corrections");
78     return false;
79   }
80   // Check that we have the vertex bias correction, needed by 
81   //   AliFMDCorrector 
82   if (what & AliForwardCorrectionManager::kVertexBias && 
83       !fcm.GetVertexBias()) { 
84     AliFatal("No event vertex bias corrections");
85     return false;
86   }
87   // Check that we have the merging efficiencies, optionally used by 
88   //   AliFMDCorrector 
89   if (what & AliForwardCorrectionManager::kMergingEfficiency && 
90       !fcm.GetMergingEfficiency()) {
91     AliFatal("No merging efficiencies");
92     return false;
93   }
94   // Check that we have the acceptance correction, needed by 
95   //   AliFMDCorrector 
96   if (what & AliForwardCorrectionManager::kAcceptance && 
97       !fcm.GetAcceptance()) { 
98     AliFatal("No acceptance corrections");
99     return false;
100   }
101   return true;
102 }
103 //____________________________________________________________________
104 Bool_t
105 AliForwardMultiplicityBase::ReadCorrections(const TAxis*& pe, 
106                                             const TAxis*& pv, 
107                                             Bool_t        mc)
108 {
109   UInt_t what = AliForwardCorrectionManager::kAll;
110   if (!fEnableLowFlux)
111     what ^= AliForwardCorrectionManager::kDoubleHit;
112   if (!GetCorrections().IsUseVertexBias())
113     what ^= AliForwardCorrectionManager::kVertexBias;
114   if (!GetCorrections().IsUseAcceptance())
115     what ^= AliForwardCorrectionManager::kAcceptance;
116   if (!GetCorrections().IsUseMergingEfficiency())
117     what ^= AliForwardCorrectionManager::kMergingEfficiency;
118
119   AliForwardCorrectionManager& fcm = AliForwardCorrectionManager::Instance();
120   if (!fcm.Init(GetEventInspector().GetCollisionSystem(),
121                 GetEventInspector().GetEnergy(),
122                 GetEventInspector().GetField(),
123                 mc,
124                 what)) return false;
125   if (!CheckCorrections(what)) return false;
126
127   // Sett our persistency pointer 
128   // fCorrManager = &fcm;
129
130   // Get the eta axis from the secondary maps - if read in
131   if (!pe) {
132     pe = fcm.GetEtaAxis();
133     if (!pe) AliFatal("No eta axis defined");
134   }
135   // Get the vertex axis from the secondary maps - if read in
136   if (!pv) {
137     pv = fcm.GetVertexAxis();
138     if (!pv) AliFatal("No vertex axis defined");
139   }
140
141   return true;
142 }
143 //____________________________________________________________________
144 AliESDEvent*
145 AliForwardMultiplicityBase::GetESDEvent()
146 {
147   AliESDEvent* esd = dynamic_cast<AliESDEvent*>(InputEvent());
148   if (!esd) {
149     AliWarning("No ESD event found for input event");
150     return 0;
151   }
152
153   // On the first event, initialize the parameters
154   if (fFirstEvent && esd->GetESDRun()) {
155     GetEventInspector().ReadRunDetails(esd);
156
157     AliInfo(Form("Initializing with parameters from the ESD:\n"
158                  "         AliESDEvent::GetBeamEnergy()   ->%f\n"
159                  "         AliESDEvent::GetBeamType()     ->%s\n"
160                  "         AliESDEvent::GetCurrentL3()    ->%f\n"
161                  "         AliESDEvent::GetMagneticField()->%f\n"
162                  "         AliESDEvent::GetRunNumber()    ->%d\n",
163                  esd->GetBeamEnergy(),
164                  esd->GetBeamType(),
165                  esd->GetCurrentL3(),
166                  esd->GetMagneticField(),
167                  esd->GetRunNumber()));
168
169     fFirstEvent = false;
170
171     InitializeSubs();
172   }
173   return esd;
174 }
175 //____________________________________________________________________
176 void
177 AliForwardMultiplicityBase::MarkEventForStore() const
178 {
179   // Make sure the AOD tree is filled 
180   AliAnalysisManager* am = AliAnalysisManager::GetAnalysisManager();
181   AliAODHandler*      ah = 
182     dynamic_cast<AliAODHandler*>(am->GetOutputEventHandler());
183   if (!ah)  
184     AliFatal("No AOD output handler set in analysis manager");
185
186   ah->SetFillAOD(kTRUE);
187 }
188
189 //____________________________________________________________________
190 void
191 AliForwardMultiplicityBase::Print(Option_t* option) const
192 {
193   // 
194   // Print information 
195   // 
196   // Parameters:
197   //    option Not used
198   //
199   
200   std::cout << "AliForwardMultiplicityBase: " << GetName() << "\n" 
201             << "  Enable low flux code:   " << (fEnableLowFlux ? "yes" : "no") 
202             << "\n"
203             << "  Off-line trigger mask:  0x" 
204             << std::hex     << std::setfill('0') 
205             << std::setw (8) << fOfflineTriggerMask 
206             << std::dec     << std::setfill (' ') << std::endl;
207   gROOT->IncreaseDirLevel();
208   if (fCorrManager) fCorrManager->Print();
209   else  
210     std::cout << "  Correction manager not set yet" << std::endl;
211   GetEventInspector()   .Print(option);
212   GetEnergyFitter()     .Print(option);    
213   GetSharingFilter()    .Print(option);
214   GetDensityCalculator().Print(option);
215   GetCorrections()      .Print(option);
216   GetHistCollector()    .Print(option);
217   gROOT->DecreaseDirLevel();
218 }
219
220
221 //
222 // EOF
223 //