]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisMuMuBase.cxx
AliAODEvent::GetHeader now return AliVHeader
[u/mrichter/AliRoot.git] / PWG / muon / AliAnalysisMuMuBase.cxx
1 #include "AliAnalysisMuMuBase.h"
2
3 /**
4  *
5  * \ingroup pwg-muon-mumu
6  *
7  * \class AliAnalysisMuMuBase
8  *
9  * Defines the interface of a sub-analysis to be steered by AliAnalysisTaskMuMu
10  *
11  * Daugther class has to implement one method :
12  *
13  * - \ref DefineHistogramCollection
14  *
15  * It may implement one or several of the FillHistosForXXX methods :
16  *
17  * - \ref FillHistosForEvent to fill histograms for one event
18  * - \ref FillHistosForMCEvent to fill histograms for one MC event
19  * - \ref FillHistosForTrack to fill histograms for one track
20  * - \ref FillHistosForPair to fill histograms for one muon track pair
21  *
22  * More rarely it may also implement :
23  *
24  * - \ref SetRun to do some changes when run number changes
25  * - \ref Terminate to finalize before ending
26  * - \ref SetEvent to e.g. append information to VEvent
27  *
28  * Daugther class can use the following methods :
29  *
30  * - \ref HasMC to know if MC information is available in the analyzed data
31  * - \ref Event to access the current event
32  * - \ref MCEvent to access to current MC event (if available)
33  *
34  * A few trivial cut methods (\ref AlwaysTrue and \ref AlwaysFalse) are defined as well and
35  * can be used to register some control cut combinations (see \ref AliAnalysisMuMuCutCombination)
36  *
37  */
38
39 #include "AliMergeableCollection.h"
40 #include "AliCounterCollection.h"
41 #include "TList.h"
42 #include "TObjString.h"
43 #include "TMath.h"
44 #include "TObjArray.h"
45 #include "AliLog.h"
46 #include "AliAnalysisMuMuBinning.h"
47 #include "TH1F.h"
48 #include "TH2F.h"
49 #include "TProfile.h"
50 #include "TRegexp.h"
51 #include "AliVEvent.h"
52 #include "AliAODEvent.h"
53 #include "AliESDEvent.h"
54 #include "AliLog.h"
55 #include "AliAnalysisMuMuCutCombination.h"
56 #include "AliAnalysisMuMuCutRegistry.h"
57
58 ClassImp(AliAnalysisMuMuBase)
59
60 namespace
61 {
62     const char* MCINPUTPREFIX = "MCINPUT";
63 };
64
65 //_____________________________________________________________________________
66 AliAnalysisMuMuBase::AliAnalysisMuMuBase()
67 :
68 TObject(),
69 fEventCounters(0x0),
70 fHistogramCollection(0x0),
71 fBinning(0x0),
72 fCutRegistry(0x0),
73 fEvent(0x0),
74 fMCEvent(0x0),
75 fHistogramToDisable(0x0),
76 fHasMC(kFALSE)
77 {
78  /// default ctor
79 }
80
81 //_____________________________________________________________________________
82 void
83 AliAnalysisMuMuBase::CreateEventHistos(UInt_t dataType,
84                                        const char* what,
85                                        const char* hname, const char* htitle,
86                                        Int_t nbinsx, Double_t xmin, Double_t xmax,
87                                        Int_t nbinsy, Double_t ymin, Double_t ymax) const
88 {
89   /** Append histograms at the event level. Depending on dataType the created histograms
90    * are for real data only, MC data only, or both
91    */
92   
93   if ( IsHistogramDisabled(hname) ) return;
94   
95   TObjArray pathNames;
96   pathNames.SetOwner(kTRUE);
97   
98   if ( dataType & kHistoForData )
99   {
100     pathNames.Add(new TObjString(Form("/%s",what)));
101   }
102   if ( ( dataType & kHistoForMCInput ) && HasMC() )
103   {
104     pathNames.Add(new TObjString(Form("/%s/%s",MCINPUTPREFIX,what)));
105   }
106   
107   CreateHistos(pathNames,hname,htitle,nbinsx,xmin,xmax,nbinsy,ymin,ymax);
108 }
109
110 //_____________________________________________________________________________
111 void
112 AliAnalysisMuMuBase::CreateEventHistos(UInt_t dataType,
113                                        const char* eventSelection,
114                                        const char* triggerClassName,
115                                        const char* centrality,
116                                        const char* hname, const char* htitle,
117                                        Int_t nbinsx, Double_t xmin, Double_t xmax,
118                                        Int_t nbinsy, Double_t ymin, Double_t ymax) const
119 {
120   /** Append histograms at the event level. Depending on dataType the created histograms
121    * are for real data only, MC data only, or both
122    */
123   
124   if ( IsHistogramDisabled(hname) ) return;
125   
126   TObjArray pathNames;
127   pathNames.SetOwner(kTRUE);
128
129   if ( dataType & kHistoForData )
130   {
131     pathNames.Add(new TObjString(Form("/%s/%s/%s",eventSelection,triggerClassName,centrality)));
132   }
133   if ( ( dataType & kHistoForMCInput ) && HasMC() )
134   {
135     pathNames.Add(new TObjString(Form("/%s/%s/%s/%s",MCINPUTPREFIX,eventSelection,triggerClassName,centrality)));
136   }
137   
138   CreateHistos(pathNames,hname,htitle,nbinsx,xmin,xmax,nbinsy,ymin,ymax);
139 }
140
141 //_____________________________________________________________________________
142 void
143 AliAnalysisMuMuBase::CreateHistos(const TObjArray& paths,
144                                   const char* hname, const char* htitle,
145                                   Int_t nbinsx, Double_t xmin, Double_t xmax,
146                                   Int_t nbinsy, Double_t ymin, Double_t ymax) const
147 {
148   /// Create multiple copies of histogram hname, one per path
149   
150   if ( IsHistogramDisabled(hname) ) return;
151
152   StdoutToAliDebug(1,paths.Print(););
153   
154   TIter next(&paths);
155   TObjString* pathName;
156
157   while ( ( pathName = static_cast<TObjString*>(next()) ) )
158   {
159     TH1* h(0x0);
160     
161     if ( nbinsy > 0 )
162     {
163       AliDebug(1,Form("Created TH2F %s/%s",pathName->String().Data(),hname));
164       h = new TH2F(hname,htitle,nbinsx,xmin,xmax,nbinsy,ymin,ymax);
165     }
166     else if ( nbinsy == 0 )
167     {
168       AliDebug(1,Form("Created TProfile %s/%s",pathName->String().Data(),hname));
169       h = new TProfile(hname,htitle,nbinsx,xmin,xmax,ymin,ymax);
170       h->Sumw2();
171     }
172     else
173     {
174       AliDebug(1,Form("Created TH1F %s/%s",pathName->String().Data(),hname));
175       h = new TH1F(hname,htitle,nbinsx,xmin,xmax);
176       
177       if ( nbinsy < -1 )
178       {
179         h->Sumw2();
180       }
181     }
182     
183     HistogramCollection()->Adopt(pathName->String().Data(),h);
184   }
185   
186   StdoutToAliDebug(1,HistogramCollection()->Print("*"););
187 }
188
189 //_____________________________________________________________________________
190 void
191 AliAnalysisMuMuBase::CreateTrackHistos(const char* eventSelection,
192                                        const char* triggerClassName,
193                                        const char* centrality,
194                                        const char* hname, const char* htitle,
195                                        Int_t nbinsx, Double_t xmin, Double_t xmax,
196                                        Int_t nbinsy, Double_t ymin, Double_t ymax) const
197 {
198   /// Create n copies of an histogram (with name=hname, title=htitle, etc..)
199   /// where n = # of track cut combinations defined
200   /// see also CreateHistos
201   
202   
203   if ( IsHistogramDisabled(hname) ) return;
204   
205   TObjArray pathNames;
206   pathNames.SetOwner(kTRUE);
207   
208   TIter nextCutCombination(CutRegistry()->GetCutCombinations(AliAnalysisMuMuCutElement::kTrack));
209   AliAnalysisMuMuCutCombination* cutCombination;
210   
211   while ( ( cutCombination = static_cast<AliAnalysisMuMuCutCombination*>(nextCutCombination())) )
212   {
213     pathNames.Add(new TObjString(Form("/%s/%s/%s/%s",eventSelection,triggerClassName,centrality,cutCombination->GetName())));
214   }
215   
216   CreateHistos(pathNames,hname,htitle,nbinsx,xmin,xmax,nbinsy,ymin,ymax);
217 }
218
219 //_____________________________________________________________________________
220 void
221 AliAnalysisMuMuBase::CreatePairHistos(const char* eventSelection,
222                                        const char* triggerClassName,
223                                        const char* centrality,
224                                        const char* hname, const char* htitle,
225                                        Int_t nbinsx, Double_t xmin, Double_t xmax,
226                                        Int_t nbinsy, Double_t ymin, Double_t ymax) const
227 {
228   /// Create n copies of an histogram (with name=hname, title=htitle, etc..)
229   /// where n = # of track *pair* cut combinations defined
230   /// see also CreateHistos
231   
232   
233   if ( IsHistogramDisabled(hname) ) return;
234   
235   TObjArray pathNames;
236   pathNames.SetOwner(kTRUE);
237   
238   TIter nextCutCombination(CutRegistry()->GetCutCombinations(AliAnalysisMuMuCutElement::kTrackPair));
239   AliAnalysisMuMuCutCombination* cutCombination;
240   
241   while ( ( cutCombination = static_cast<AliAnalysisMuMuCutCombination*>(nextCutCombination())) )
242   {
243     pathNames.Add(new TObjString(Form("/%s/%s/%s/%s",eventSelection,triggerClassName,centrality,cutCombination->GetName())));
244   }
245   
246   CreateHistos(pathNames,hname,htitle,nbinsx,xmin,xmax,nbinsy,ymin,ymax);
247 }
248
249 //_____________________________________________________________________________
250 void AliAnalysisMuMuBase::DisableHistograms(const char* pattern)
251 {
252   /// Disable the histogramming of all the histograms matching the pattern
253   
254   TString spattern(pattern);
255   if (spattern=="*")
256   {
257     delete fHistogramToDisable;
258     fHistogramToDisable = 0x0;
259   }
260   
261   if (!fHistogramToDisable)
262   {
263     fHistogramToDisable = new TList;
264     fHistogramToDisable->SetOwner(kTRUE);
265   }
266   
267   fHistogramToDisable->Add(new TObjString(spattern));
268 }
269
270 //_____________________________________________________________________________
271 Int_t AliAnalysisMuMuBase::GetNbins(Double_t xmin, Double_t xmax, Double_t xstep)
272 {
273   /// Compute number of bins to get a given step between each values in the xmin,xmax range
274   if ( TMath::AreEqualRel(xstep,0.0,1E-9) ) return 1;
275   
276   return TMath::Nint(TMath::Abs((xmax-xmin)/xstep));
277 }
278
279 //_____________________________________________________________________________
280 TH1* AliAnalysisMuMuBase::Histo(const char* eventSelection, const char* triggerClassName, const char* histoname)
281 {
282   /// Get one histo back
283   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s",eventSelection,triggerClassName,histoname)) : 0x0;
284 }
285
286 //_____________________________________________________________________________
287 TH1* AliAnalysisMuMuBase::Histo(const char* eventSelection, const char* histoname)
288 {
289   /// Get one histo back
290   return fHistogramCollection ? fHistogramCollection->Histo(eventSelection,histoname) : 0x0;
291 }
292
293 //_____________________________________________________________________________
294 TH1* AliAnalysisMuMuBase::Histo(const char* eventSelection,
295                                 const char* triggerClassName,
296                                 const char* cent,
297                                 const char* histoname)
298 {
299   /// Get one histo back
300   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s",eventSelection,triggerClassName,cent),histoname) : 0x0;
301 }
302
303 //_____________________________________________________________________________
304 TH1* AliAnalysisMuMuBase::Histo(const char* eventSelection,
305                                 const char* triggerClassName,
306                                 const char* cent,
307                                 const char* what,
308                                 const char* histoname)
309 {
310   /// Get one histo back
311   
312   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s/%s",eventSelection,triggerClassName,cent,what),histoname) : 0x0;
313 }
314
315 //_____________________________________________________________________________
316 void AliAnalysisMuMuBase::Init(AliCounterCollection& cc,
317                                AliMergeableCollection& hc,
318                                const AliAnalysisMuMuBinning& binning,
319                                const AliAnalysisMuMuCutRegistry& registry)
320 {
321   /// Set the internal references
322   fEventCounters = &cc;
323   fHistogramCollection = &hc;
324   fBinning = &binning;
325   fCutRegistry = &registry;
326 }
327
328 //_____________________________________________________________________________
329 Bool_t AliAnalysisMuMuBase::IsHistogramDisabled(const char* hname) const
330 {
331   /// Whether or not a given histogram (identified by its name)
332   /// is disabled or not
333   
334   if ( !fHistogramToDisable )
335   {
336     return kFALSE;
337   }
338   TString shname(hname);
339   TIter next(fHistogramToDisable);
340   TObjString* str(0x0);
341   
342   while ( ( str = static_cast<TObjString*>(next()) ) )
343   {
344     if ( shname.Contains(TRegexp(str->String()) ) )
345     {
346       return kTRUE;
347     }
348   }
349   return kFALSE;
350 }
351
352 //_____________________________________________________________________________
353 Bool_t AliAnalysisMuMuBase::IsHistogrammingDisabled() const
354 {
355   /// Whether or not *all* histograms are disabled
356   
357   if ( fHistogramToDisable && fHistogramToDisable->GetEntries()==1 )
358   {
359     TObjString* r = static_cast<TObjString*>(fHistogramToDisable->First());
360     if ( r->String() == "*" )
361     {
362       return kTRUE;
363     }
364   }
365   return kFALSE;
366 }
367
368 //_____________________________________________________________________________
369 TH1* AliAnalysisMuMuBase::MCHisto(const char* eventSelection, const char* triggerClassName, const char* histoname)
370 {
371   /// Get one histo back
372   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s/%s",MCINPUTPREFIX,eventSelection,triggerClassName,histoname)) : 0x0;
373 }
374
375 //_____________________________________________________________________________
376 TH1* AliAnalysisMuMuBase::MCHisto(const char* eventSelection, const char* histoname)
377 {
378   /// Get one histo back
379   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s",MCINPUTPREFIX,eventSelection,histoname)) : 0x0;
380 }
381
382 //_____________________________________________________________________________
383 TH1* AliAnalysisMuMuBase::MCHisto(const char* eventSelection,
384                                   const char* triggerClassName,
385                                   const char* cent,
386                                   const char* histoname)
387 {
388   /// Get one histo back
389   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s/%s",MCINPUTPREFIX,eventSelection,triggerClassName,cent),histoname) : 0x0;
390 }
391
392 //_____________________________________________________________________________
393 TH1* AliAnalysisMuMuBase::MCHisto(const char* eventSelection,
394                                   const char* triggerClassName,
395                                   const char* cent,
396                                   const char* what,
397                                   const char* histoname)
398 {
399   /// Get one histo back
400   
401   return fHistogramCollection ? fHistogramCollection->Histo(Form("/%s/%s/%s/%s/%s",MCINPUTPREFIX,eventSelection,triggerClassName,cent,what),histoname) : 0x0;
402 }
403
404 //_____________________________________________________________________________
405 void AliAnalysisMuMuBase::SetEvent(AliVEvent* event, AliMCEvent* mcEvent)
406 {
407   /// Set pointers to event (and mcEvent)
408   /// This is the place where a sub analysis might (if really needed)
409   /// append things to the event. Please DO NOT alter the original content,
410   /// only add to it. Unless of course you want to shoot yourself in the
411   /// foot...
412   
413   fEvent = event;
414   fMCEvent = mcEvent;
415 }
416