]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisMuMuEventCutter.cxx
Split the TaskMuMu into more manageable sub-analysis (Laurent)
[u/mrichter/AliRoot.git] / PWG / muon / AliAnalysisMuMuEventCutter.cxx
1 #include "AliAnalysisMuMuEventCutter.h"
2
3 /**
4  *
5  * \ingroup pwg-muon-mumu
6  *
7  * \class AliAnalysisMuMuEventCutter
8  *
9  * We're grouping here various event cut methods that can be used
10  * as AliAnalysisMuMuCutElement. For instance :
11  *
12  * - \ref IsPhysicsSelected this is the normal physics selection check
13  * - \ref IsPhysicsSelectedVDM version of the physics selection used for VdM analysis
14  * - \ref IsSPDzVertexInRange selects event with a SPDvertex in a given range
15  * - \ref IsAbsZBelowValue selects event with |Zvertex| (being SPD or not) below a given value
16  * - \ref IsAbsZBelowValue as above but requesting explicitely the SPD vertex
17  * - \ref IsSPDzQA whether the SPD vertex is a good one
18  */
19
20 #include "TObjString.h"
21 #include "AliLog.h"
22 #include "AliMuonEventCuts.h"
23 #include "TList.h"
24 #include "Riostream.h"
25 #include "AliVVertex.h"
26 #include "AliAODVertex.h"
27 #include "AliVVZERO.h"
28 #include "AliInputEventHandler.h"
29 #include "AliVEvent.h"
30 #include "TMath.h"
31 #include "AliESDTZERO.h"
32 #include "AliAODTZERO.h"
33 #include "AliESDEvent.h"
34 #include "AliAODEvent.h"
35
36 ClassImp(AliAnalysisMuMuEventCutter)
37
38 //______________________________________________________________________________
39 AliAnalysisMuMuEventCutter::AliAnalysisMuMuEventCutter(TList* triggerClasses)
40 : TObject(), fMuonEventCuts(0x0)
41 {
42   /// ctor
43   TString tclasses;
44   
45   if ( !triggerClasses )
46   {
47     tclasses = "ANY";
48   }
49   else
50   {
51     TObjString* tname;
52     TIter next(triggerClasses);
53     
54     while ( ( tname = static_cast<TObjString*>(next()) ) )
55     {
56       if (tclasses.Length()>0)
57       {
58         tclasses += ",";
59       }
60       
61       tclasses += tname->String();
62     }
63   }
64   
65   MuonEventCuts()->SetTrigClassPatterns(tclasses);
66 }
67
68 //______________________________________________________________________________
69 AliAnalysisMuMuEventCutter::~AliAnalysisMuMuEventCutter()
70 {
71   /// dtor
72   delete fMuonEventCuts;
73 }
74
75 //_____________________________________________________________________________
76 Bool_t AliAnalysisMuMuEventCutter::SelectTriggerClass(const TString& firedTriggerClasses,
77                                                       TString& acceptedClasses,
78                                                       UInt_t L0, UInt_t L1, UInt_t L2) const
79 {
80   /// Forward the trigger class selection to MuonEventCuts::GetSelectedTrigClassesInEvent
81   acceptedClasses = "";
82   
83   TIter next(MuonEventCuts()->GetSelectedTrigClassesInEvent(firedTriggerClasses,L0,L1,L2));
84   TObjString* str;
85   
86   while ( ( str = static_cast<TObjString*>(next()) ) )
87   {
88     acceptedClasses += str->String();
89     acceptedClasses += " ";
90   }
91   return (acceptedClasses.Length()>0);
92 }
93
94 //_____________________________________________________________________________
95 Bool_t AliAnalysisMuMuEventCutter::IsPhysicsSelected(const AliInputEventHandler& eventHandler) const
96 {
97   /// Whether or not the event is physics selected
98   return const_cast<AliInputEventHandler&>(eventHandler).IsEventSelected() & AliVEvent::kAny;
99 }
100
101 //_____________________________________________________________________________
102 Bool_t AliAnalysisMuMuEventCutter::IsPhysicsSelectedVDM(const AliVEvent& event) const
103 {
104   // cut used in vdM scans
105   
106   AliVVZERO* vzero = event.GetVZEROData();
107   
108   if (vzero)
109   {
110     Float_t v0a = vzero->GetV0ATime();
111     Float_t v0c = vzero->GetV0CTime();
112     
113     Float_t v0diff = v0a-v0c;
114     Float_t v0sum = v0a+v0c;
115     
116     if ( ( v0sum > 10.5 && v0sum < 18 ) && ( v0diff > 4 && v0diff < 12 ) )
117     {
118       return kTRUE;
119     }
120   }
121   return kFALSE;
122 }
123
124 //_____________________________________________________________________________
125 Bool_t AliAnalysisMuMuEventCutter::IsAbsZBelowValue(const AliVEvent& event, const Double_t& z) const
126 {
127   // Checks if the absolute value of the Z component of the primary vertex is below a certain value
128   
129   const AliVVertex* vertex = event.GetPrimaryVertex();
130   return (TMath::Abs(vertex->GetZ())<=z);
131 }
132
133 //_____________________________________________________________________________
134 Bool_t AliAnalysisMuMuEventCutter::IsAbsZSPDBelowValue(const AliVEvent& event, const Double_t& z) const
135 {
136   // Checks if the absolute value of the SPD Z component of the given vertex is below a certain value
137   AliAODVertex* SPDVertex = static_cast<const AliAODEvent*>(&event)->GetPrimaryVertexSPD();
138   
139   if ( !SPDVertex )
140   {
141     AliError("SPD |z| cut requested and no SPD vertex found in the event");
142     return kFALSE;
143   }
144   
145   else return (TMath::Abs(SPDVertex->GetZ())<=z);
146 }
147
148
149 //_____________________________________________________________________________
150 Bool_t AliAnalysisMuMuEventCutter::IsSPDzVertexInRange(AliVEvent& event, const Double_t& zMin, const Double_t& zMax) const
151 {
152   /// Whether or not the SPD Z vertex is in the range [zMin,zMax[
153   
154   AliAODVertex* SPDVertex = static_cast<AliAODEvent*>(&event)->GetPrimaryVertexSPD();
155   
156   if ( !SPDVertex )
157   {
158     AliError("Cut on SPD z Vertex requested for an event with no SPD vertex info");
159     return kFALSE;
160   }
161   Double_t zV = SPDVertex->GetZ();
162   
163   if ( zV >= zMin && zV < zMax ) return kTRUE;
164   else return kFALSE;
165 }
166
167
168 //_____________________________________________________________________________
169 Bool_t AliAnalysisMuMuEventCutter::HasSPDVertex(AliVEvent& event) const
170 {
171   /// Does the event have a SPD vertex ?
172   AliAODVertex* SPDVertex = static_cast<AliAODEvent*>(&event)->GetPrimaryVertexSPD();
173   if ( SPDVertex && SPDVertex->GetNContributors() > 0) return kTRUE;
174   else return kFALSE;
175 }
176
177
178 //_____________________________________________________________________________
179 Bool_t AliAnalysisMuMuEventCutter::IsSPDzQA(const AliVEvent& event, /*const AliVVertex& vertex2Test,*/ const Double_t& zResCut, const Double_t& zDifCut) const
180 {
181   // Checks if the value of the Z component of the given vertex fullfills the quality assurance condition
182   
183   const AliVVertex* vertex = event.GetPrimaryVertex();
184   const AliAODVertex* vertex2Test = static_cast<const AliAODEvent*>(&event)->GetPrimaryVertexSPD();
185   
186   Double_t cov[6]={0};
187   vertex2Test->GetCovarianceMatrix(cov);
188   Double_t zRes = TMath::Sqrt(cov[5]);
189   
190   if ( (zRes < zResCut) && TMath::Abs(vertex2Test->GetZ() - vertex->GetZ()) <= zDifCut )
191   {
192     return kTRUE;
193   }
194   else return kFALSE;
195 }
196
197
198 //_____________________________________________________________________________
199 Bool_t AliAnalysisMuMuEventCutter::IsMeandNchdEtaInRange(AliVEvent& event, const Double_t& dNchdEtaMin, const Double_t& dNchdEtaMax) const
200 {
201   TList* nchList = static_cast<TList*>(event.FindListObject("NCH"));
202   
203   if (!nchList || nchList->IsEmpty())
204   {
205     AliFatal("No NCH information found in event. Nch analysis MUST be executed to apply a NCH cut");
206     return kFALSE;
207   }
208   
209   Int_t i(0);
210   
211   while ( nchList->At(i)->IsA() != TObjString::Class() ) //Asign a name to find it by name
212   {
213     i++;
214   }
215   
216   TObjString* value = static_cast<TObjString*>(nchList->At(i));
217   Double_t meandNchdEta = value->String().Atof();
218   
219   if ( meandNchdEta >= dNchdEtaMin && meandNchdEta < dNchdEtaMax ) return kTRUE;
220   else return kFALSE;
221 }
222
223
224 //_____________________________________________________________________________
225 AliMuonEventCuts*
226 AliAnalysisMuMuEventCutter::MuonEventCuts() const
227 {
228   /// Return the single instance of AliMuonEventCuts object we're using
229   
230   if (!fMuonEventCuts)
231   {
232     fMuonEventCuts = new AliMuonEventCuts("EventCut","");
233   }
234   return fMuonEventCuts;
235 }
236
237
238 //_____________________________________________________________________________
239 void AliAnalysisMuMuEventCutter::NameOfIsSPDzVertexInRange(TString& name, const Double_t& zMin, const Double_t& zMax) const
240 {
241   name.Form("SPDZBTW%3.2f_%3.2f",zMin,zMax);
242 }
243 //_____________________________________________________________________________
244 void AliAnalysisMuMuEventCutter::NameOfIsAbsZBelowValue(TString& name, const Double_t& z) const
245 {
246   name.Form("ABSZLT%3.2f",z);
247 }
248
249 //_____________________________________________________________________________
250 void AliAnalysisMuMuEventCutter::NameOfIsAbsZSPDBelowValue(TString& name, const Double_t& z) const
251 {
252   name.Form("SPDABSZLT%3.2f",z);
253 }
254
255
256 //_____________________________________________________________________________
257 void AliAnalysisMuMuEventCutter::NameOfIsSPDzQA(TString& name, const Double_t& zResCut, const Double_t& zDifCut) const
258 {
259   name.Form("SPDZQA_RES%3.2f_ZDIF%3.2f",zResCut,zDifCut);
260 }
261
262 //_____________________________________________________________________________
263 void AliAnalysisMuMuEventCutter::NameOfIsMeandNchdEtaInRange(TString& name, const Double_t& dNchdEtaMin, const Double_t& dNchdEtaMax) const
264 {
265   name.Form("MEANDNDETABTW%3.2f_%3.2f",dNchdEtaMin,dNchdEtaMax);
266 }
267
268
269 //_____________________________________________________________________________
270 Bool_t AliAnalysisMuMuEventCutter::IsTZEROPileUp(const AliVEvent& event) const
271 {
272   Bool_t pileupFlag(kFALSE);
273   
274   if ( event.IsA() == AliESDEvent::Class() )
275   {
276     const AliESDTZERO* tzero = static_cast<AliESDEvent&>(const_cast<AliVEvent&>(event)).GetESDTZERO();
277     if ( tzero )
278     {
279       pileupFlag = tzero->GetPileupFlag();
280     }
281   }
282   else if ( event.IsA() == AliAODEvent::Class() )
283   {
284     AliAODTZERO* tzero = static_cast<const AliAODEvent&>(event).GetTZEROData();
285     if ( tzero )
286     {
287       pileupFlag = tzero->GetPileupFlag();
288     }
289   }
290   return pileupFlag;
291 }