]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisTaskSingleMu.cxx
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWG / muon / AliAnalysisTaskSingleMu.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 //-----------------------------------------------------------------------------
19 /// \class AliAnalysisTaskSingleMu
20 /// Analysis task for single muons in the spectrometer.
21 /// The output is a list of histograms and CF containers.
22 /// The macro class can run on AODs or ESDs.
23 /// If Monte Carlo information is present, some basics checks are performed.
24 ///
25 /// \author Diego Stocco
26 //-----------------------------------------------------------------------------
27
28 #define AliAnalysisTaskSingleMu_cxx
29
30 #include "AliAnalysisTaskSingleMu.h"
31
32 // ROOT includes
33 #include "TROOT.h"
34 #include "TH1.h"
35 #include "TH2.h"
36 #include "TAxis.h"
37 #include "TCanvas.h"
38 #include "TLegend.h"
39 #include "TMath.h"
40 #include "TObjString.h"
41 #include "TObjArray.h"
42 #include "TF1.h"
43 #include "TStyle.h"
44 //#include "TMCProcess.h"
45 #include "TArrayI.h"
46 #include "TPaveStats.h"
47 #include "TFitResultPtr.h"
48 #include "TFile.h"
49 #include "THashList.h"
50
51 // STEER includes
52 #include "AliAODEvent.h"
53 #include "AliAODTrack.h"
54 #include "AliAODMCParticle.h"
55 #include "AliMCEvent.h"
56 #include "AliMCParticle.h"
57 #include "AliESDEvent.h"
58 #include "AliESDMuonTrack.h"
59 #include "AliVHeader.h"
60 #include "AliAODMCHeader.h"
61 #include "AliStack.h"
62
63 // ANALYSIS includes
64 #include "AliAnalysisManager.h"
65
66 // CORRFW includes
67 #include "AliCFContainer.h"
68 #include "AliCFGridSparse.h"
69 #include "AliCFEffGrid.h"
70
71 // PWG includes
72 #include "AliVAnalysisMuon.h"
73 #include "AliMergeableCollection.h"
74 #include "AliCounterCollection.h"
75 #include "AliMuonEventCuts.h"
76 #include "AliMuonTrackCuts.h"
77 #include "AliAnalysisMuonUtility.h"
78
79
80 /// \cond CLASSIMP
81 ClassImp(AliAnalysisTaskSingleMu) // Class implementation in ROOT context
82 /// \endcond
83
84
85 //________________________________________________________________________
86 AliAnalysisTaskSingleMu::AliAnalysisTaskSingleMu() :
87   AliVAnalysisMuon(),
88   fThetaAbsKeys(0x0),
89   fCutOnDimu(kFALSE)
90 {
91   /// Default ctor.
92 }
93
94 //________________________________________________________________________
95 AliAnalysisTaskSingleMu::AliAnalysisTaskSingleMu(const char *name, const AliMuonTrackCuts& cuts) :
96   AliVAnalysisMuon(name, cuts),
97   fThetaAbsKeys(0x0),
98   fCutOnDimu(kFALSE)
99 {
100   //
101   /// Constructor.
102   //
103   TString thetaAbsKeys = "ThetaAbs23 ThetaAbs310";
104   fThetaAbsKeys = thetaAbsKeys.Tokenize(" ");
105 }
106
107
108 //________________________________________________________________________
109 AliAnalysisTaskSingleMu::~AliAnalysisTaskSingleMu()
110 {
111   //
112   /// Destructor
113   //
114
115   delete fThetaAbsKeys;
116 }
117
118 //___________________________________________________________________________
119 void AliAnalysisTaskSingleMu::MyUserCreateOutputObjects()
120 {
121
122   TH1* histo = 0x0;
123   TString histoName = "", histoTitle = "";
124   
125   Int_t nVzBins = 40;
126   Double_t vzMin = -20., vzMax = 20.;
127   TString vzName("Vz"), vzTitle("Vz"), vzUnits("cm");  
128   
129   histoName = "hIpVtx";
130   histo = new TH1F(histoName.Data(), histoName.Data(), nVzBins, vzMin, vzMax);
131   histo->SetXTitle("v_{z} (cm)");
132   AddObjectToCollection(histo, kIPVz);
133
134   
135   Int_t nPtBins = 80;
136   Double_t ptMin = 0., ptMax = 80.;
137   TString ptName("Pt"), ptTitle("p_{T}"), ptUnits("GeV/c");
138   
139   Int_t nEtaBins = 25;
140   Double_t etaMin = -4.5, etaMax = -2.;
141   TString etaName("Eta"), etaTitle("#eta"), etaUnits("");
142   
143   Int_t nPhiBins = 36;
144   Double_t phiMin = 0.; Double_t phiMax = 2.*TMath::Pi();
145   TString phiName("Phi"), phiTitle("#phi"), phiUnits("rad");
146     
147   Int_t nChargeBins = 2;
148   Double_t chargeMin = -2, chargeMax = 2.;
149   TString chargeName("Charge"), chargeTitle("charge"), chargeUnits("e");
150   
151   Int_t nThetaAbsEndBins = 2;
152   Double_t thetaAbsEndMin = -0.5, thetaAbsEndMax = 1.5;
153   TString thetaAbsEndName("ThetaAbsEnd"), thetaAbsEndTitle("#theta_{abs}"), thetaAbsEndUnits("a.u.");    
154   
155   Int_t nMotherTypeBins = kNtrackSources;
156   Double_t motherTypeMin = -0.5, motherTypeMax = (Double_t)kNtrackSources - 0.5;
157   TString motherType("MotherType"), motherTypeTitle("motherType"), motherTypeUnits("");
158     
159   Int_t nbins[kNvars] = {nPtBins, nEtaBins, nPhiBins, nVzBins, nChargeBins, nThetaAbsEndBins, nMotherTypeBins};
160   Double_t xmin[kNvars] = {ptMin, etaMin, phiMin, vzMin, chargeMin, thetaAbsEndMin, motherTypeMin};
161   Double_t xmax[kNvars] = {ptMax, etaMax, phiMax, vzMax, chargeMax, thetaAbsEndMax, motherTypeMax};
162   TString axisTitle[kNvars] = {ptTitle, etaTitle, phiTitle, vzTitle, chargeTitle, thetaAbsEndTitle, motherTypeTitle};
163   TString axisUnits[kNvars] = {ptUnits, etaUnits, phiUnits, vzUnits, chargeUnits, thetaAbsEndUnits, motherTypeUnits};
164
165   AliCFContainer* cfContainer = new AliCFContainer("SingleMuContainer","Container for tracks",kNsteps,kNvars,nbins);
166   
167   for ( Int_t idim = 0; idim<kNvars; idim++){
168     histoTitle = Form("%s (%s)", axisTitle[idim].Data(), axisUnits[idim].Data());
169     histoTitle.ReplaceAll("()","");
170     
171     cfContainer->SetVarTitle(idim, histoTitle.Data());
172     cfContainer->SetBinLimits(idim, xmin[idim], xmax[idim]);
173   }
174   
175   TString stepTitle[kNsteps] = {"reconstructed", "generated"};
176
177   TAxis* currAxis = 0x0;
178   for (Int_t istep=0; istep<kNsteps; istep++){
179     cfContainer->SetStepTitle(istep, stepTitle[istep].Data());
180     AliCFGridSparse* gridSparse = cfContainer->GetGrid(istep);
181         
182     currAxis = gridSparse->GetAxis(kHvarMotherType);
183     for ( Int_t ibin=0; ibin<fSrcKeys->GetEntries(); ibin++ ) {
184       currAxis->SetBinLabel(ibin+1, fSrcKeys->At(ibin)->GetName());
185     }
186   }
187   
188   AddObjectToCollection(cfContainer, kTrackContainer);
189   
190   histoName = "hRecoDimu";
191   TH2* histoDimu = new TH2F(histoName.Data(), histoName.Data(), 24, 0., 120., 30, 0., 150.);
192   histoDimu->SetXTitle("min. p_{T}^{#mu} (GeV/c)");
193   histoDimu->SetYTitle("M_{#mu#mu} (GeV/c^{2})");
194   AddObjectToCollection(histoDimu, kNobjectTypes);
195   
196   histoName = "hGeneratedZ";
197   AddObjectToCollection(histoDimu->Clone(histoName.Data()), kNobjectTypes+1);
198   
199   
200   histoName = "hZmuEtaCorr";
201   histoDimu = new TH2F(histoName.Data(), histoName.Data(), 160, -8., 8., 160,-8., 8.);
202   histoDimu->SetXTitle("#eta_{#mu-}");
203   histoDimu->SetYTitle("#eta_{#mu+}");
204   AddObjectToCollection(histoDimu, kNobjectTypes+2);
205   
206   fMuonTrackCuts->Print("mask");
207 }
208
209 //________________________________________________________________________
210 void AliAnalysisTaskSingleMu::ProcessEvent(TString physSel, const TObjArray& selectTrigClasses, TString centrality)
211 {
212   //
213   /// Fill output objects
214   //
215
216   Double_t ipVz = AliAnalysisMuonUtility::GetVertexSPD(InputEvent())->GetZ();
217   Double_t ipVzMC = MCEvent() ? AliAnalysisMuonUtility::GetMCVertexZ(InputEvent(),MCEvent()) : 0.;
218   
219   for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
220     TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
221     ((TH1*)GetMergeableObject(physSel, trigClassName, centrality, "hIpVtx"))->Fill(ipVz);
222   }
223
224   // Bool_t isPileupFromSPD = ( fAODEvent && ! fAODEvent->GetTracklets() ) ? InputEvent()->IsPileupFromSPD(3, 0.8, 3., 2., 5.) : InputEvent()->IsPileupFromSPDInMultBins(); // Avoid break when reading Muon AODs (tracklet info is not present and IsPileupFromSPDInMultBins crashes // UNCOMMENT TO REJECT PILEUP
225   // if ( isPileupFromSPD ) return; // UNCOMMENT TO REJECT PILEUP
226   
227   Double_t containerInput[kNvars];
228   AliVParticle* track = 0x0;
229
230   Int_t nSteps = MCEvent() ? 2 : 1;
231   for ( Int_t istep = 0; istep<nSteps; ++istep ) {
232     Int_t nTracks = ( istep == kStepReconstructed ) ? AliAnalysisMuonUtility::GetNTracks(InputEvent()) : MCEvent()->GetNumberOfTracks();
233     
234     TObjArray selectedTracks(nTracks);
235     TArrayI trackSources(nTracks);
236     for (Int_t itrack = 0; itrack < nTracks; itrack++) {
237       track = ( istep == kStepReconstructed ) ? AliAnalysisMuonUtility::GetTrack(itrack,InputEvent()) : MCEvent()->GetTrack(itrack);
238       
239       // In case of MC we usually ask that the particle is a muon
240       // However, in W or Z simulations, Pythia stores both the initial muon
241       // (before ISR, FSR and kt kick) and the final state one.
242       // The first muon is of course there only for information and should be rejected.
243       // The Pythia code for initial state particles is 21
244       Bool_t isSelected = ( istep == kStepReconstructed ) ? fMuonTrackCuts->IsSelected(track) : ( TMath::Abs(track->PdgCode()) == 13 && AliAnalysisMuonUtility::GetStatusCode(track) != 21 );
245       if ( ! isSelected ) continue;
246       
247       selectedTracks.AddAt(track,itrack);
248       trackSources[itrack] = GetParticleType(track);
249       
250     } // loop on tracks
251     
252     // Loop on selected tracks
253     TArrayI rejectTrack(nTracks);
254     for ( Int_t itrack=0; itrack<nTracks; itrack++) {
255       track = static_cast<AliVParticle*>(selectedTracks.At(itrack));
256       if ( ! track ) continue;
257       
258       if ( istep == kStepGeneratedMC || fCutOnDimu ) {
259         // Check dimuons
260         for ( Int_t jtrack=itrack+1; jtrack<nTracks; jtrack++ ) {
261           AliVParticle* auxTrack = static_cast<AliVParticle*>(selectedTracks.At(jtrack));
262           if ( ! auxTrack ) continue;
263           if ( track->Charge() * auxTrack->Charge() >= 0 ) continue;
264           
265           TLorentzVector dimuPair = AliAnalysisMuonUtility::GetTrackPair(track,auxTrack);
266           Double_t ptMin = TMath::Min(track->Pt(),auxTrack->Pt());
267           Double_t invMass = dimuPair.M();
268           if ( istep == kStepReconstructed ) {
269             if ( invMass > 60. && invMass < 120. ) {
270               rejectTrack[itrack] = 1;
271               rejectTrack[jtrack] = 1;
272             }
273             for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
274               TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
275               if ( ! fMuonTrackCuts->TrackPtCutMatchTrigClass(track, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) || ! fMuonTrackCuts->TrackPtCutMatchTrigClass(auxTrack, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) ) continue;
276               ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hRecoDimu"))->Fill(ptMin,invMass);
277             } // loop on triggers
278           }
279           else {
280             if ( trackSources[itrack] == kZbosonMu && trackSources[jtrack] == kZbosonMu ) {
281               Bool_t isAccepted = kTRUE;
282               AliVParticle* muDaughter[2] = {0x0, 0x0};
283               for ( Int_t imu=0; imu<2; imu++ ) {
284                 AliVParticle* currPart = ( imu == 0 ) ? track : auxTrack;
285                 if ( currPart->Charge() < 0. ) muDaughter[0] = currPart;
286                 else muDaughter[1] = currPart;
287                 if ( currPart->Eta() < -4.5 || currPart->Eta() > -2. ) {
288                   isAccepted = kFALSE;
289                 }
290               } // loop on muons in the pair
291               
292               Double_t pairRapidity = dimuPair.Rapidity();
293               if ( pairRapidity < -4. || pairRapidity > -2.5 ) isAccepted = kFALSE;
294               //printf("Rapidity Z %g  pair %g\n",track->Y(), pairRapidity);
295               
296               for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
297                 TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
298                 if ( isAccepted )  ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hGeneratedZ"))->Fill(ptMin,invMass);
299                 ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hZmuEtaCorr"))->Fill(muDaughter[0]->Eta(),muDaughter[1]->Eta());
300               } // loop on selected trig
301             } // both muons from Z
302           } // kStepGeneratedMC
303         } // loop on auxiliary tracks
304       } // apply cut on dimu
305       if ( rejectTrack[itrack] > 0 ) continue;
306       
307       Double_t thetaAbsEndDeg = 0;
308       if ( istep == kStepReconstructed ) {
309         thetaAbsEndDeg = AliAnalysisMuonUtility::GetThetaAbsDeg(track);
310       }
311       else {
312         thetaAbsEndDeg = ( TMath::Pi()-track->Theta() ) * TMath::RadToDeg();
313       }
314       Int_t thetaAbsBin = ( thetaAbsEndDeg < 3. ) ? kThetaAbs23 : kThetaAbs310;
315
316       containerInput[kHvarPt]         = track->Pt();
317       containerInput[kHvarEta]        = track->Eta();
318       containerInput[kHvarPhi]        = track->Phi();
319       containerInput[kHvarVz]         = ( istep == kStepReconstructed ) ? ipVz : ipVzMC;
320       containerInput[kHvarCharge]     = track->Charge()/3.;
321       containerInput[kHvarThetaAbs]   = (Double_t)thetaAbsBin;
322       containerInput[kHvarMotherType] = (Double_t)trackSources[itrack];
323       
324       for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
325         TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
326         if ( istep == kStepReconstructed && ! fMuonTrackCuts->TrackPtCutMatchTrigClass(track, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) ) continue;
327         ((AliCFContainer*)GetMergeableObject(physSel, trigClassName, centrality, "SingleMuContainer"))->Fill(containerInput,istep);
328       } // loop on selected trigger classes
329     } // loop on tracks
330   } // loop on container steps
331 }
332
333
334 //________________________________________________________________________
335 void AliAnalysisTaskSingleMu::Terminate(Option_t *) {
336   //
337   /// Draw some histograms at the end.
338   //
339   
340   AliVAnalysisMuon::Terminate("");
341   
342   if ( ! fMergeableCollection ) return;
343   
344   TString physSel = fTerminateOptions->At(0)->GetName();
345   TString trigClassName = fTerminateOptions->At(1)->GetName();
346   TString centralityRange = fTerminateOptions->At(2)->GetName();
347   TString furtherOpt = fTerminateOptions->At(3)->GetName();
348   
349   TString minBiasTrig = "";
350   TObjArray* optArr = furtherOpt.Tokenize(" ");
351   TString currName = "";
352   for ( Int_t iopt=0; iopt<optArr->GetEntries(); iopt++ ) {
353     currName = optArr->At(iopt)->GetName();
354     if ( currName.Contains("-B-") || currName.Contains("ANY") ) minBiasTrig = currName;
355   }
356   delete optArr;
357   
358   furtherOpt.ToUpper();
359   Bool_t plotChargeAsymmetry = furtherOpt.Contains("ASYM");
360   
361   AliCFContainer* cfContainer = static_cast<AliCFContainer*> ( GetSum(physSel,trigClassName,centralityRange,"SingleMuContainer") );
362   if ( ! cfContainer ) return;
363   
364   AliCFEffGrid* effSparse = new AliCFEffGrid(Form("eff%s", cfContainer->GetName()),Form("Efficiency %s", cfContainer->GetTitle()),*cfContainer);
365   effSparse->CalculateEfficiency(kStepReconstructed, kStepGeneratedMC);
366   
367   AliCFGridSparse* gridSparseArray[3] = {effSparse->GetNum(), effSparse->GetDen(), effSparse};
368   TString gridSparseName[3] = {cfContainer->GetStepTitle(kStepReconstructed), cfContainer->GetStepTitle(kStepGeneratedMC), "Efficiency"};
369   
370   Int_t srcColors[kNtrackSources] = {kBlack, kRed, kSpring, kTeal, kBlue, kViolet, kMagenta, kOrange, kGray};
371   //  TString allSrcNames = "";
372   //  for ( Int_t isrc=0; isrc<kNtrackSources; ++isrc ) {
373   //    if ( ! allSrcNames.IsNull() ) allSrcNames.Append(" ");
374   //    allSrcNames += fSrcKeys->At(isrc)->GetName();
375   //  }
376   
377   TCanvas* can = 0x0;
378   Int_t xshift = 100;
379   Int_t yshift = 100;
380   Int_t igroup1 = -1;
381   Int_t igroup2 = 0;
382   
383   Bool_t isMC = furtherOpt.Contains("MC");
384   
385   TAxis* srcAxis = gridSparseArray[0]->GetAxis(kHvarMotherType);
386   Int_t unIdBin = srcAxis->GetNbins();
387   for ( Int_t ibin=1; ibin<=srcAxis->GetNbins(); ibin++ ) {
388     currName = srcAxis->GetBinLabel(ibin);
389     if ( currName.Contains("Unidentified") ) unIdBin = ibin;
390   }
391   
392   Int_t firstSrcBin = ( isMC ) ? 1 : unIdBin;
393   Int_t lastSrcBin  = ( isMC ) ? srcAxis->GetNbins() - 1 : unIdBin;
394   if ( ! isMC ) srcColors[unIdBin-1] = 1;
395   
396   ////////////////
397   // Kinematics //
398   ////////////////
399   TString chargeNames[3] = {fChargeKeys->At(0)->GetName(), fChargeKeys->At(1)->GetName(), "Total"};
400   THashList histoList[3];
401   for ( Int_t icharge=0; icharge<3; icharge++ ) {
402     histoList[icharge].SetName(chargeNames[icharge].Data());
403   }
404   for ( Int_t isrc = firstSrcBin; isrc <= lastSrcBin; ++isrc ) {
405     for ( Int_t icharge=0; icharge<3; ++icharge ) {
406       Int_t icharge1 = ( icharge < 2 ) ? icharge : 0;
407       Int_t icharge2 = ( icharge < 2 ) ? icharge : 1;
408       for ( Int_t igrid=0; igrid<3; ++igrid ) {
409         if ( gridSparseArray[igrid]->GetEntries() == 0. ) break;
410         if ( gridSparseArray[igrid]->IsA() != AliCFEffGrid::Class() ) {
411           SetSparseRange(gridSparseArray[igrid], kHvarEta, "", -3.999, -2.501);
412           SetSparseRange(gridSparseArray[igrid], kHvarMotherType, "", isrc, isrc, "USEBIN");
413           SetSparseRange(gridSparseArray[igrid], kHvarCharge, "", icharge1+1, icharge2+1, "USEBIN");
414         }
415         TH1* histo = gridSparseArray[igrid]->Project(kHvarPt, kHvarEta);
416         histo->SetName(Form("hPtEta_%s_%s_%s", gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), chargeNames[icharge].Data()));
417         if ( histo->Integral() > 0 ) histoList[icharge].Add(histo);
418         for ( Int_t iproj=0; iproj<4; ++iproj ) {
419           histo = gridSparseArray[igrid]->Project(iproj);
420           histo->SetName(Form("proj%i_%s_%s_%s", iproj, gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), chargeNames[icharge].Data()));
421           if ( histo->Integral() > 0 ) histoList[icharge].Add(histo);
422         } // loop on projections
423       } // loop on grid sparse
424     } // loop on charge
425   } // loop on track sources
426   
427   // Get charge asymmetry or mu+/mu-
428   THashList histoListRatio;
429   TString basePlotName = plotChargeAsymmetry ? "ChargeAsymmetry" : "ChargeRatio";
430   histoListRatio.SetName(basePlotName.Data());
431   Int_t baseCharge = 1;
432   Int_t auxCharge = 1-baseCharge;
433   for ( Int_t ihisto=0; ihisto<histoList[baseCharge].GetEntries(); ihisto++ ) {
434     TObject* obj = histoList[baseCharge].At(ihisto);
435     TString histoName = obj->GetName();
436     if ( histoName.Contains(gridSparseName[2].Data()) ) continue;
437     TString searchName = histoName;
438     searchName.ReplaceAll(fChargeKeys->At(baseCharge)->GetName(), fChargeKeys->At(auxCharge)->GetName());
439     TH1* auxHisto = static_cast<TH1*> (histoList[auxCharge].FindObject(searchName.Data()));
440     if ( ! auxHisto ) continue;
441     histoName.ReplaceAll(fChargeKeys->At(baseCharge)->GetName(),basePlotName.Data());
442     TH1* histo = static_cast<TH1*> (obj->Clone(histoName.Data()));
443     if ( plotChargeAsymmetry ) {
444       histo->Add(auxHisto, -1.);
445       // h2 + h1 = 2xh2 + (h1-h2)
446       auxHisto->Add(auxHisto, histo, 2.);
447     }
448     histo->Divide(auxHisto);
449     TString axisTitle = plotChargeAsymmetry ? Form("(%s - %s) / (%s + %s)", fChargeKeys->At(1)->GetName(), fChargeKeys->At(0)->GetName(), fChargeKeys->At(1)->GetName(), fChargeKeys->At(0)->GetName()) : Form("%s / %s", fChargeKeys->At(1)->GetName(), fChargeKeys->At(0)->GetName());
450     axisTitle.ReplaceAll("MuPlus","#mu^{+}");
451     axisTitle.ReplaceAll("MuMinus","#mu^{-}");
452     histo->GetYaxis()->SetTitle(axisTitle.Data());
453     histo->SetStats(kFALSE);
454     histoListRatio.Add(histo);
455   }
456   
457   // Plot kinematics
458   TString histoName = "", drawOpt = "";
459   for ( Int_t itype=0; itype<3; itype++ ) {
460     THashList* currList = 0x0;
461     Int_t nCharges = 1;
462     if ( itype == 1 ) currList = &histoListRatio;
463     else if ( itype == 2 ) currList = &histoList[2];
464     else nCharges = 2;
465     for ( Int_t igrid=0; igrid<3; ++igrid ) {
466       igroup1 = igrid;
467       TCanvas* canKine = 0x0;
468       TLegend* legKine = 0x0;
469       for ( Int_t iproj=0; iproj<4; ++iproj ) {
470         for ( Int_t isrc = firstSrcBin; isrc <= lastSrcBin; ++isrc ) {
471           for ( Int_t icharge=0; icharge<nCharges; ++icharge ) {
472             if ( itype == 0 ) currList = &histoList[icharge];
473             histoName = Form("proj%i_%s_%s_%s", iproj, gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), currList->GetName());
474             TH1* histo = static_cast<TH1*>(currList->FindObject(histoName.Data()));
475             if ( ! histo ) continue;
476             if ( ! canKine ) {
477               igroup2 = itype;
478               igroup1 = igrid;
479               currName = Form("%s_%s_%s", GetName(), currList->GetName(), gridSparseName[igrid].Data());
480               canKine = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
481               canKine->Divide(2,2);
482               legKine = new TLegend(0.6, 0.6, 0.8, 0.8);
483             }
484             canKine->cd(iproj+1);
485             if ( itype != 1 ) {
486               if ( ( iproj == kHvarPt || iproj == kHvarVz ) && gridSparseArray[igrid]->IsA() != AliCFEffGrid::Class() ) gPad->SetLogy();
487             }
488             Bool_t isFirst = ( gPad->GetListOfPrimitives()->GetEntries() == 0 );
489             drawOpt = isFirst ? "e" : "esames";
490             histo->SetLineColor(srcColors[isrc-1]);
491             histo->SetMarkerColor(srcColors[isrc-1]);
492             histo->SetMarkerStyle(20+4*icharge);
493             histo->Draw(drawOpt.Data());
494             TPaveStats* paveStats = (TPaveStats*)histo->FindObject("stats");
495             if ( paveStats ) paveStats->SetTextColor(srcColors[isrc-1]);
496             if ( iproj == 0 ) {
497               TString legEntry = ( itype == 0 ) ? fChargeKeys->At(icharge)->GetName() : "";
498               if ( isMC ) legEntry += Form(" %s", srcAxis->GetBinLabel(isrc));
499               if ( ! legEntry.IsNull() ) legKine->AddEntry(histo,legEntry.Data(), "lp");
500             }
501           } // loop on mu charge
502         } // loop on track sources
503       } // loop on projections
504       if ( legKine && legKine->GetNRows() > 0 ) {
505         canKine->cd(1);
506         legKine->Draw("same");
507       }
508     } // loop on grid sparse
509   } // loop on types
510   
511   
512   for ( Int_t igrid=0; igrid<3; igrid++ ) {
513     if ( gridSparseArray[igrid]->IsA() == AliCFEffGrid::Class() ) continue;
514     SetSparseRange(gridSparseArray[igrid], kHvarCharge, "", 1, gridSparseArray[igrid]->GetAxis(kHvarCharge)->GetNbins(), "USEBIN"); // Reset range
515   } // loop on container steps
516   
517   //////////////////////
518   // Event statistics //
519   //////////////////////
520   printf("\nTotal analyzed events:\n");
521   TString evtSel = Form("trigger:%s", trigClassName.Data());
522   fEventCounters->PrintSum(evtSel.Data());
523   printf("Physics selected analyzed events:\n");
524   evtSel = Form("trigger:%s/selected:yes", trigClassName.Data());
525   fEventCounters->PrintSum(evtSel.Data());
526   
527   TString countPhysSel = "any";
528   if ( physSel.Contains(fPhysSelKeys->At(kPhysSelPass)->GetName()) ) countPhysSel = "yes";
529   else if ( physSel.Contains(fPhysSelKeys->At(kPhysSelReject)->GetName()) ) countPhysSel="no";
530   countPhysSel.Prepend("selected:");
531   printf("Analyzed events vs. centrality:\n");
532   evtSel = Form("trigger:%s/%s", trigClassName.Data(), countPhysSel.Data());
533   fEventCounters->Print("centrality",evtSel.Data(),kTRUE);
534   
535   
536   TString outFilename = Form("/tmp/out%s.root", GetName());
537   TFile* outFile = new TFile(outFilename.Data(),"RECREATE");
538   for ( Int_t icharge=0; icharge<3; icharge++ ) {
539     histoList[icharge].Write();
540   }
541   histoListRatio.Write();
542   outFile->Close();
543   printf("\nCreating file %s\n", outFilename.Data(
544          ));
545   
546   ///////////////////
547   // Vertex method //
548   ///////////////////
549   if ( ! furtherOpt.Contains("VERTEX") ) return;
550   Int_t firstMother = ( isMC ) ? 0 : kUnidentified;
551   Int_t lastMother = ( isMC ) ? kNtrackSources - 1 : kUnidentified;
552   igroup1++;
553   TH1* eventVertex = (TH1*)GetSum(physSel, minBiasTrig, centralityRange, "hIpVtx");
554   if ( ! eventVertex ) return;
555   Double_t minZ = -9.99, maxZ = 9.99;
556   Double_t meanZ = 0., sigmaZ = 4.;
557   Double_t nSigma = 2.;
558   TString fitOpt = "R0S";
559   Bool_t fixFitRange = kFALSE;
560   TString fitFormula = Form("[0]+[1]*(x+[2])");
561   
562   // Get vertex shape
563   if ( eventVertex->GetSumw2N() == 0 ) eventVertex->Sumw2();
564   Double_t eventVtxIntegral = eventVertex->Integral(0,eventVertex->GetNbinsX()+1); // Include under/overflow
565   printf("Event vertex integral %.0f\n\n", eventVtxIntegral);
566   if ( eventVtxIntegral <= 0. ) return;
567   eventVertex->Scale(1./eventVtxIntegral);
568   printf("\nFit MB vertex\n");
569   eventVertex->Fit("gaus",fitOpt.Data(),"",minZ,maxZ);
570   TF1* vtxFit = (TF1*)eventVertex->GetListOfFunctions()->FindObject("gaus");
571   currName = "vtxIntegrated";
572   can = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
573   can->SetLogy();
574   eventVertex->Draw();
575   vtxFit->Draw("same");
576   
577   
578   enum {kRecoHF, kRecoBkg, kInputHF, kInputDecay, kRecoAll, kNrecoHistos};
579   TString baseRecoName[kNrecoHistos] = {"RecoHF", "RecoBkg", "InputHF", "InputDecay", "RecoAll"};
580   TArrayI sumMothers[kNrecoHistos];
581   sumMothers[kRecoHF].Set(0);
582   sumMothers[kRecoBkg].Set(0);
583   sumMothers[kInputHF].Set(3);
584   sumMothers[kInputHF][0] = kCharmMu;
585   sumMothers[kInputHF][1] = kBeautyMu;
586   sumMothers[kInputHF][2] = kQuarkoniumMu;
587   sumMothers[kInputDecay].Set(1);
588   sumMothers[kInputDecay][0] = kDecayMu;
589   sumMothers[kRecoAll].Set(srcAxis->GetNbins());
590   for ( Int_t isrc=1; isrc<srcAxis->GetNbins(); ++isrc ) {
591     sumMothers[kRecoAll][isrc-1] = isrc;
592   }
593   
594   meanZ = vtxFit->GetParameter(1);
595   sigmaZ = vtxFit->GetParameter(2);
596   
597   Double_t minZfit = ( fixFitRange ) ? minZ : meanZ - nSigma*sigmaZ;
598   Double_t maxZfit = ( fixFitRange ) ? maxZ : meanZ + nSigma*sigmaZ;
599   
600   TF1* fitFunc = new TF1("fitFunc", fitFormula.Data(), minZ, maxZ);
601   fitFunc->SetLineColor(2);
602   fitFunc->SetParNames("Line norm", "Line slope", "Free path");
603   const Double_t kFreePath = 153.; // 150.; // 130.; // cm
604   //fitFunc->SetParameters(0.,1.);
605   fitFunc->FixParameter(2, kFreePath);
606   
607   AliCFGridSparse* gridSparse = cfContainer->GetGrid(kStepReconstructed);
608   TAxis* ptAxis = gridSparse->GetAxis(kHvarPt);
609   
610   Double_t slope = 0.;
611   Double_t limitNorm = 0., limitSlope = 0.;
612   Int_t firstPtBin = 0, lastPtBin = 0;
613   
614   gStyle->SetOptFit(1111);
615   
616   for ( Int_t itheta=0; itheta<kNthetaAbs; ++itheta ) {
617     igroup2++;
618     SetSparseRange(gridSparse, kHvarThetaAbs, "", itheta+1, itheta+1, "USEBIN");
619     SetSparseRange(gridSparse, kHvarPt, "", 1, ptAxis->GetNbins(), "USEBIN");
620     TH1* recoHisto[kNrecoHistos];
621     for ( Int_t ireco=0; ireco<kNrecoHistos; ++ireco ) {
622       recoHisto[ireco] = gridSparse->Project(kHvarPt);
623       histoName = Form("%sMuon_%s", baseRecoName[ireco].Data(), fThetaAbsKeys->At(itheta)->GetName());
624       recoHisto[ireco]->SetName(histoName.Data());
625       recoHisto[ireco]->SetTitle(histoName.Data());
626       recoHisto[ireco]->Reset();
627       recoHisto[ireco]->Sumw2();
628       for ( Int_t isrc=0; isrc<sumMothers[ireco].GetSize(); ++isrc ) {
629         SetSparseRange(gridSparse, kHvarMotherType, "", sumMothers[ireco][isrc], sumMothers[ireco][isrc], "USEBIN");
630         TH1* auxHisto = gridSparse->Project(kHvarPt);
631         recoHisto[ireco]->Add(auxHisto);
632         delete auxHisto;
633       }
634     }
635     SetSparseRange(gridSparse, kHvarMotherType, "", firstMother+1, lastMother+1, "USEBIN");
636     Int_t currDraw = 0;
637     
638     for ( Int_t ibinpt=0; ibinpt<=ptAxis->GetNbins(); ++ibinpt ) {
639       firstPtBin = ibinpt;
640       lastPtBin = ( ibinpt == 0 ) ? ptAxis->GetNbins() : ibinpt;
641       SetSparseRange(gridSparse, kHvarPt, "", firstPtBin, lastPtBin, "USEBIN");
642       TH1* histo = gridSparse->Project(kHvarVz);
643       histo->SetName(Form("hVtx_%s_%s_ptBin%i", cfContainer->GetStepTitle(kStepReconstructed), fThetaAbsKeys->At(itheta)->GetName(), ibinpt));
644       if ( histo->Integral() < 100. ) break;
645       printf("\nFit %.2f < pt < %.2f (entries %g)\n", ptAxis->GetBinLowEdge(firstPtBin), ptAxis->GetBinUpEdge(lastPtBin), histo->GetEntries());
646       histo->Divide(eventVertex);
647       Double_t norm = histo->GetBinContent(histo->FindBin(0.));
648       histo->GetYaxis()->SetTitle("#frac{dN_{#mu}}{dv_{z}} / #left(#frac{1}{N_{MB}}#frac{dN_{MB}}{dv_{z}}#right)");
649       slope = ( histo->GetBinContent(histo->FindBin(meanZ+sigmaZ)) -
650                histo->GetBinContent(histo->FindBin(meanZ-sigmaZ)) ) / ( 2. * sigmaZ );
651       
652       if ( slope < 0. ) slope = norm/kFreePath;
653       
654       // Try to fit twice: it fit fails the first time
655       // set some limits on parameters
656       for ( Int_t itry=0; itry<2; itry++ ) {
657         fitFunc->SetParameter(0, norm);
658         fitFunc->SetParameter(1, slope);
659         if ( itry > 0 ) {
660           limitNorm = 2.*histo->Integral();
661           limitSlope = 2.*histo->Integral()/kFreePath;
662           //fitFunc->SetParLimits(0, 0., limitNorm); // REMEMBER TO CHECK
663           fitFunc->SetParLimits(1, 0., limitSlope); // REMEMBER TO CHECK
664           printf("Norm 0. < %f < %f  slope  0. < %f < %f\n", norm, limitNorm, slope, limitSlope);
665         }
666         TFitResultPtr fitRes = histo->Fit(fitFunc, fitOpt.Data(), "", minZfit, maxZfit);
667         
668         //      if ( gMinuit->fCstatu.Contains("CONVERGED") &&
669         if ( ((Int_t)fitRes) == 0 &&
670             fitFunc->GetParameter(0) > 0. &&
671             fitFunc->GetParameter(1) > 0. )
672           break;
673         else if ( furtherOpt.Contains("REFIT") ) printf("Re-fit with limits\n");
674         else {
675           printf("Warning: fit problems !!!\n");
676           break;
677         }
678       }
679       
680       Double_t p0 = fitFunc->GetParameter(0);
681       Double_t p0err = fitFunc->GetParError(0);
682       Double_t p1 = fitFunc->GetParameter(1);
683       Double_t p1err = fitFunc->GetParError(1);
684       
685       Double_t nullVz = ( p1 != 0. ) ? -p0/p1 : 0.;
686       Double_t nullVzErr = ( p0 != 0. && p1 != 0. ) ? TMath::Abs(nullVz) * TMath::Sqrt(p0err*p0err/(p0*p0) + p1err*p1err/(p1*p1) ) : 0.;
687       
688       printf("Null value at %f +- %f\n", nullVz - kFreePath, nullVzErr);
689       
690       recoHisto[kRecoHF]->SetBinContent(ibinpt, p0);
691       recoHisto[kRecoHF]->SetBinError(ibinpt, p0err);
692       recoHisto[kRecoBkg]->SetBinContent(ibinpt, ( kFreePath + meanZ ) * p1);
693       recoHisto[kRecoBkg]->SetBinError(ibinpt, ( kFreePath + meanZ ) * p1err);
694       if ( currDraw%4 == 0 ){
695         currName = Form("vtx_%s_PtBin%i",fThetaAbsKeys->At(itheta)->GetName(), ibinpt);
696         can = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
697         can->Divide(2,2);
698       }
699       can->cd( currDraw%4 + 1 );
700       can->SetLogy();
701       histo->Draw();
702       fitFunc->DrawCopy("same");
703       currDraw++;
704     } // loop on pt bins
705     SetSparseRange(gridSparse, kHvarMotherType, "", firstMother+1, lastMother+1, "USEBIN");
706     currName = Form("recoPt_%s",fThetaAbsKeys->At(itheta)->GetName());
707     can = new TCanvas(currName.Data(),currName.Data(),(igroup1+1)*xshift,igroup2*yshift,600,600);
708     TLegend* leg = new TLegend(0.6, 0.6, 0.8, 0.8);
709     drawOpt = "e";
710     for ( Int_t ireco=0; ireco<kNrecoHistos-1; ++ireco ) {
711       if ( recoHisto[ireco]->GetEntries() == 0. ) continue;
712       TH1* ratio = (TH1*)recoHisto[ireco]->Clone(Form("%s_ratio", recoHisto[ireco]->GetName()));
713       ratio->Divide(recoHisto[kRecoAll]);
714       ratio->SetLineColor(srcColors[ireco]);
715       ratio->SetMarkerColor(srcColors[ireco]);
716       ratio->SetMarkerStyle(20+ireco);
717       ratio->GetYaxis()->SetTitle("fraction of total");
718       ratio->Draw(drawOpt.Data());
719       leg->AddEntry(ratio,baseRecoName[ireco].Data(), "lp");
720       drawOpt = "esame";
721     }
722     leg->Draw("same");
723   } // loop on theta abs
724 }