]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/muon/AliAnalysisTaskSingleMu.cxx
Fix treatment of MC sources in terminate (with Vertex method) (Diego)
[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   TH2* histoGenZ = static_cast<TH2*>(histoDimu->Clone(histoName.Data()));
198   histoGenZ->SetTitle(histoName.Data());
199   AddObjectToCollection(histoGenZ, kNobjectTypes+1);
200   
201   
202   histoName = "hZmuEtaCorr";
203   histoDimu = new TH2F(histoName.Data(), histoName.Data(), 160, -8., 8., 160, -8., 8.);
204   histoDimu->SetXTitle("#eta_{#mu-}");
205   histoDimu->SetYTitle("#eta_{#mu+}");
206   AddObjectToCollection(histoDimu, kNobjectTypes+2);
207   
208   fMuonTrackCuts->Print("mask");
209   
210   AliInfo(Form("Apply cut on dimuon (60<M_mumu<120 GeV/c^2) to reject Z contribution: %i", fCutOnDimu));
211 }
212
213 //________________________________________________________________________
214 void AliAnalysisTaskSingleMu::ProcessEvent(TString physSel, const TObjArray& selectTrigClasses, TString centrality)
215 {
216   //
217   /// Fill output objects
218   //
219
220   Double_t ipVz = AliAnalysisMuonUtility::GetVertexSPD(InputEvent())->GetZ();
221   Double_t ipVzMC = MCEvent() ? AliAnalysisMuonUtility::GetMCVertexZ(InputEvent(),MCEvent()) : 0.;
222   
223   for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
224     TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
225     ((TH1*)GetMergeableObject(physSel, trigClassName, centrality, "hIpVtx"))->Fill(ipVz);
226   }
227   
228   TF1* beautyMuWgt = static_cast<TF1*>(GetWeight("beautyMu"));
229
230   // 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
231   // if ( isPileupFromSPD ) return; // UNCOMMENT TO REJECT PILEUP
232   
233   Double_t containerInput[kNvars];
234   AliVParticle* track = 0x0;
235
236   Int_t nSteps = MCEvent() ? 2 : 1;
237   for ( Int_t istep = 0; istep<nSteps; ++istep ) {
238     Int_t nTracks = ( istep == kStepReconstructed ) ? AliAnalysisMuonUtility::GetNTracks(InputEvent()) : MCEvent()->GetNumberOfTracks();
239     
240     TObjArray selectedTracks(nTracks);
241     TArrayI trackSources(nTracks);
242     TArrayF trackWgt(nTracks);
243     trackWgt.Reset(1.);
244     for (Int_t itrack = 0; itrack < nTracks; itrack++) {
245       track = ( istep == kStepReconstructed ) ? AliAnalysisMuonUtility::GetTrack(itrack,InputEvent()) : MCEvent()->GetTrack(itrack);
246       
247       // In case of MC we usually ask that the particle is a muon
248       // However, in W or Z simulations, Pythia stores both the initial muon
249       // (before ISR, FSR and kt kick) and the final state one.
250       // The first muon is of course there only for information and should be rejected.
251       // The Pythia code for initial state particles is 21
252       Bool_t isSelected = ( istep == kStepReconstructed ) ? fMuonTrackCuts->IsSelected(track) : ( TMath::Abs(track->PdgCode()) == 13 && AliAnalysisMuonUtility::GetStatusCode(track) != 21 );
253       if ( ! isSelected ) continue;
254       
255       selectedTracks.AddAt(track,itrack);
256       trackSources[itrack] = GetParticleType(track);
257       
258       if ( trackSources[itrack] == kBeautyMu && beautyMuWgt ) {
259         AliVParticle* mcTrack = ( kStepReconstructed ) ? MCEvent()->GetTrack(track->GetLabel()) : track;
260         trackWgt[itrack] = beautyMuWgt->Eval(mcTrack->Pt());
261 //        Int_t iancestor = fUtilityMuonAncestor->GetAncestor(track,MCEvent());
262 //        AliVParticle* motherPart = MCEvent()->GetTrack(iancestor);
263 //        trackWgt[itrack] = beautyMuWgt->GetBinContent(beautyMuWgt->GetXaxis()->FindBin(motherPart->Pt()));
264       }
265     } // loop on tracks
266     
267     // Loop on selected tracks
268     TArrayI rejectTrack(nTracks);
269     for ( Int_t itrack=0; itrack<nTracks; itrack++) {
270       track = static_cast<AliVParticle*>(selectedTracks.At(itrack));
271       if ( ! track ) continue;
272       
273       // Check dimuons
274       for ( Int_t jtrack=itrack+1; jtrack<nTracks; jtrack++ ) {
275         AliVParticle* auxTrack = static_cast<AliVParticle*>(selectedTracks.At(jtrack));
276         if ( ! auxTrack ) continue;
277         if ( track->Charge() * auxTrack->Charge() >= 0 ) continue;
278           
279         TLorentzVector dimuPair = AliAnalysisMuonUtility::GetTrackPair(track,auxTrack);
280         Double_t ptMin = TMath::Min(track->Pt(),auxTrack->Pt());
281         Double_t invMass = dimuPair.M();
282         if ( fCutOnDimu && invMass > 60. && invMass < 120. ) {
283           rejectTrack[itrack] = 1;
284           rejectTrack[jtrack] = 1;
285         }
286         Double_t dimuWgt = trackWgt[itrack] * trackWgt[jtrack];
287         if ( istep == kStepReconstructed ) {
288           for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
289             TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
290             if ( ! fMuonTrackCuts->TrackPtCutMatchTrigClass(track, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) || ! fMuonTrackCuts->TrackPtCutMatchTrigClass(auxTrack, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) ) continue;
291             ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hRecoDimu"))->Fill(ptMin,invMass,dimuWgt);
292           } // loop on triggers
293         }
294         else {
295           if ( trackSources[itrack] == kZbosonMu && trackSources[jtrack] == kZbosonMu ) {
296             Bool_t isAccepted = kTRUE;
297             AliVParticle* muDaughter[2] = {0x0, 0x0};
298             for ( Int_t imu=0; imu<2; imu++ ) {
299               AliVParticle* currPart = ( imu == 0 ) ? track : auxTrack;
300               if ( currPart->Charge() < 0. ) muDaughter[0] = currPart;
301               else muDaughter[1] = currPart;
302               if ( currPart->Eta() < -4.5 || currPart->Eta() > -2. ) {
303                 isAccepted = kFALSE;
304               }
305             } // loop on muons in the pair
306               
307             Double_t pairRapidity = dimuPair.Rapidity();
308             if ( pairRapidity < -4. || pairRapidity > -2.5 ) isAccepted = kFALSE;
309             //printf("Rapidity Z %g  pair %g\n",track->Y(), pairRapidity);
310               
311             for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
312               TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
313               if ( isAccepted ) ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hGeneratedZ"))->Fill(ptMin,invMass,dimuWgt);
314               ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hZmuEtaCorr"))->Fill(muDaughter[0]->Eta(),muDaughter[1]->Eta(),dimuWgt);
315             } // loop on selected trig
316           } // both muons from Z
317         } // kStepGeneratedMC
318       } // loop on auxiliary tracks
319       if ( rejectTrack[itrack] > 0 ) continue;
320       
321       Double_t thetaAbsEndDeg = 0;
322       if ( istep == kStepReconstructed ) {
323         thetaAbsEndDeg = AliAnalysisMuonUtility::GetThetaAbsDeg(track);
324       }
325       else {
326         thetaAbsEndDeg = ( TMath::Pi()-track->Theta() ) * TMath::RadToDeg();
327       }
328       Int_t thetaAbsBin = ( thetaAbsEndDeg < 3. ) ? kThetaAbs23 : kThetaAbs310;
329
330       containerInput[kHvarPt]         = track->Pt();
331       containerInput[kHvarEta]        = track->Eta();
332       containerInput[kHvarPhi]        = track->Phi();
333       containerInput[kHvarVz]         = ( istep == kStepReconstructed ) ? ipVz : ipVzMC;
334       containerInput[kHvarCharge]     = track->Charge()/3.;
335       containerInput[kHvarThetaAbs]   = (Double_t)thetaAbsBin;
336       containerInput[kHvarMotherType] = (Double_t)trackSources[itrack];
337       
338       for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
339         TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
340         if ( istep == kStepReconstructed && ! fMuonTrackCuts->TrackPtCutMatchTrigClass(track, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) ) continue;
341         ((AliCFContainer*)GetMergeableObject(physSel, trigClassName, centrality, "SingleMuContainer"))->Fill(containerInput,istep,trackWgt[itrack]);
342       } // loop on selected trigger classes
343     } // loop on tracks
344   } // loop on container steps
345 }
346
347
348 //________________________________________________________________________
349 void AliAnalysisTaskSingleMu::Terminate(Option_t *) {
350   //
351   /// Draw some histograms at the end.
352   //
353   
354   AliVAnalysisMuon::Terminate("");
355   
356   if ( ! fMergeableCollection ) return;
357   
358   TString physSel = fTerminateOptions->At(0)->GetName();
359   TString trigClassName = fTerminateOptions->At(1)->GetName();
360   TString centralityRange = fTerminateOptions->At(2)->GetName();
361   TString furtherOpt = fTerminateOptions->At(3)->GetName();
362   
363   TString minBiasTrig = "";
364   TObjArray* optArr = furtherOpt.Tokenize(" ");
365   TString currName = "";
366   for ( Int_t iopt=0; iopt<optArr->GetEntries(); iopt++ ) {
367     currName = optArr->At(iopt)->GetName();
368     if ( currName.Contains("-B-") || currName.Contains("ANY") ) minBiasTrig = currName;
369   }
370   delete optArr;
371   
372   furtherOpt.ToUpper();
373   Bool_t plotChargeAsymmetry = furtherOpt.Contains("ASYM");
374   
375   AliCFContainer* inContainer = static_cast<AliCFContainer*> ( GetSum(physSel,trigClassName,centralityRange,"SingleMuContainer") );
376   if ( ! inContainer ) return;
377   
378   AliCFContainer* cfContainer = inContainer;
379   
380   if ( ! furtherOpt.Contains("GENINTRIGCLASS") && trigClassName != "ANY" ) {
381     // The output container contains both the reconstructed and (in MC)
382     // the generated muons in a specific trigger class.
383     // Since the trigger pt level of the track is matched to the trigger class,
384     // analyzing the MUHigh trigger (for example) is equivalent of analyzing
385     // Hpt tracks.
386     // However, in this way, the generated muons distribution depend
387     // on a condition on the reconstructed track.
388     // If we calculate the Acc.xEff. as the ratio of reconstructed and
389     // generated muons IN A TRIGGER CLASS, we are biasing the final result.
390     // The correct value of Acc.xEff. is hence obtained as the distribution
391     // of reconstructed tracks in a muon trigger class, divided by the
392     // total number of generated class (which is in the "class" ANY).
393     // The following code sets the generated muons as the one in the class ANY.
394     // The feature is the default. If you want the generated muons in the same
395     // trigger class as the generated tracks, use option "GENINTRIGCLASS"
396     AliCFContainer* fullMCcontainer = static_cast<AliCFContainer*> ( GetSum(physSel,"ANY",centralityRange,"SingleMuContainer") );
397     if ( fullMCcontainer ) {
398       cfContainer = static_cast<AliCFContainer*>(cfContainer->Clone("SingleMuContainerCombo"));
399       cfContainer->SetGrid(kStepGeneratedMC,fullMCcontainer->GetGrid(kStepGeneratedMC));
400     }
401   }
402   
403   AliCFEffGrid* effSparse = new AliCFEffGrid(Form("eff%s", cfContainer->GetName()),Form("Efficiency %s", cfContainer->GetTitle()),*cfContainer);
404   effSparse->CalculateEfficiency(kStepReconstructed, kStepGeneratedMC);
405   
406   AliCFGridSparse* gridSparseArray[3] = {effSparse->GetNum(), effSparse->GetDen(), effSparse};
407   TString gridSparseName[3] = {cfContainer->GetStepTitle(kStepReconstructed), cfContainer->GetStepTitle(kStepGeneratedMC), "Efficiency"};
408   
409   Int_t srcColors[kNtrackSources] = {kBlack, kRed, kSpring, kTeal, kBlue, kViolet, kMagenta, kOrange, kGray};
410   //  TString allSrcNames = "";
411   //  for ( Int_t isrc=0; isrc<kNtrackSources; ++isrc ) {
412   //    if ( ! allSrcNames.IsNull() ) allSrcNames.Append(" ");
413   //    allSrcNames += fSrcKeys->At(isrc)->GetName();
414   //  }
415   
416   TCanvas* can = 0x0;
417   Int_t xshift = 100;
418   Int_t yshift = 100;
419   Int_t igroup1 = -1;
420   Int_t igroup2 = 0;
421   
422   Bool_t isMC = furtherOpt.Contains("MC");
423   
424   TAxis* srcAxis = gridSparseArray[0]->GetAxis(kHvarMotherType);
425   Int_t unIdBin = srcAxis->FindBin(fSrcKeys->At(kUnidentified)->GetName());
426   if ( unIdBin < 1 ) unIdBin = srcAxis->GetNbins();
427   
428   Int_t firstSrcBin = ( isMC ) ? 1 : unIdBin;
429   Int_t lastSrcBin  = ( isMC ) ? srcAxis->GetNbins() - 1 : unIdBin;
430   if ( ! isMC ) srcColors[unIdBin-1] = 1;
431   
432   ////////////////
433   // Kinematics //
434   ////////////////
435   TString chargeNames[3] = {fChargeKeys->At(0)->GetName(), fChargeKeys->At(1)->GetName(), "Total"};
436   THashList histoList[3];
437   for ( Int_t icharge=0; icharge<3; icharge++ ) {
438     histoList[icharge].SetName(chargeNames[icharge].Data());
439   }
440   for ( Int_t isrc = firstSrcBin; isrc <= lastSrcBin; ++isrc ) {
441     for ( Int_t icharge=0; icharge<3; ++icharge ) {
442       Int_t icharge1 = ( icharge < 2 ) ? icharge : 0;
443       Int_t icharge2 = ( icharge < 2 ) ? icharge : 1;
444       for ( Int_t igrid=0; igrid<3; ++igrid ) {
445         if ( gridSparseArray[igrid]->GetEntries() == 0. ) break;
446         if ( gridSparseArray[igrid]->IsA() != AliCFEffGrid::Class() ) {
447           SetSparseRange(gridSparseArray[igrid], kHvarEta, "", -3.999, -2.501);
448           SetSparseRange(gridSparseArray[igrid], kHvarMotherType, "", isrc, isrc, "USEBIN");
449           SetSparseRange(gridSparseArray[igrid], kHvarCharge, "", icharge1+1, icharge2+1, "USEBIN");
450         }
451         TH1* histo = gridSparseArray[igrid]->Project(kHvarPt, kHvarEta);
452         histo->SetName(Form("hPtEta_%s_%s_%s", gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), chargeNames[icharge].Data()));
453         if ( histo->Integral() > 0 ) histoList[icharge].Add(histo);
454         for ( Int_t iproj=0; iproj<4; ++iproj ) {
455           histo = gridSparseArray[igrid]->Project(iproj);
456           histo->SetName(Form("proj%i_%s_%s_%s", iproj, gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), chargeNames[icharge].Data()));
457           if ( histo->Integral() > 0 ) histoList[icharge].Add(histo);
458         } // loop on projections
459       } // loop on grid sparse
460     } // loop on charge
461   } // loop on track sources
462   
463   // Get charge asymmetry or mu+/mu-
464   THashList histoListRatio;
465   TString basePlotName = plotChargeAsymmetry ? "ChargeAsymmetry" : "ChargeRatio";
466   histoListRatio.SetName(basePlotName.Data());
467   Int_t baseCharge = 1;
468   Int_t auxCharge = 1-baseCharge;
469   for ( Int_t ihisto=0; ihisto<histoList[baseCharge].GetEntries(); ihisto++ ) {
470     TObject* obj = histoList[baseCharge].At(ihisto);
471     TString histoName = obj->GetName();
472     if ( histoName.Contains(gridSparseName[2].Data()) ) continue;
473     TString searchName = histoName;
474     searchName.ReplaceAll(fChargeKeys->At(baseCharge)->GetName(), fChargeKeys->At(auxCharge)->GetName());
475     TH1* auxHisto = static_cast<TH1*> (histoList[auxCharge].FindObject(searchName.Data()));
476     if ( ! auxHisto ) continue;
477     histoName.ReplaceAll(fChargeKeys->At(baseCharge)->GetName(),basePlotName.Data());
478     TH1* histo = static_cast<TH1*> (obj->Clone(histoName.Data()));
479     if ( plotChargeAsymmetry ) {
480       histo->Add(auxHisto, -1.);
481       // h2 + h1 = 2xh2 + (h1-h2)
482       auxHisto->Add(auxHisto, histo, 2.);
483     }
484     histo->Divide(auxHisto);
485     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());
486     axisTitle.ReplaceAll("MuPlus","#mu^{+}");
487     axisTitle.ReplaceAll("MuMinus","#mu^{-}");
488     histo->GetYaxis()->SetTitle(axisTitle.Data());
489     histo->SetStats(kFALSE);
490     histoListRatio.Add(histo);
491   }
492   
493   // Plot kinematics
494   TString histoName = "", drawOpt = "";
495   for ( Int_t itype=0; itype<3; itype++ ) {
496     THashList* currList = 0x0;
497     Int_t nCharges = 1;
498     if ( itype == 1 ) currList = &histoListRatio;
499     else if ( itype == 2 ) currList = &histoList[2];
500     else nCharges = 2;
501     for ( Int_t igrid=0; igrid<3; ++igrid ) {
502       igroup1 = igrid;
503       TCanvas* canKine = 0x0;
504       TLegend* legKine = 0x0;
505       for ( Int_t iproj=0; iproj<4; ++iproj ) {
506         for ( Int_t isrc = firstSrcBin; isrc <= lastSrcBin; ++isrc ) {
507           for ( Int_t icharge=0; icharge<nCharges; ++icharge ) {
508             if ( itype == 0 ) currList = &histoList[icharge];
509             histoName = Form("proj%i_%s_%s_%s", iproj, gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), currList->GetName());
510             TH1* histo = static_cast<TH1*>(currList->FindObject(histoName.Data()));
511             if ( ! histo ) continue;
512             if ( ! canKine ) {
513               igroup2 = itype;
514               igroup1 = igrid;
515               currName = Form("%s_%s_%s", GetName(), currList->GetName(), gridSparseName[igrid].Data());
516               canKine = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
517               canKine->Divide(2,2);
518               legKine = new TLegend(0.6, 0.6, 0.8, 0.8);
519             }
520             canKine->cd(iproj+1);
521             if ( itype != 1 ) {
522               if ( ( iproj == kHvarPt || iproj == kHvarVz ) && gridSparseArray[igrid]->IsA() != AliCFEffGrid::Class() ) gPad->SetLogy();
523             }
524             Bool_t isFirst = ( gPad->GetListOfPrimitives()->GetEntries() == 0 );
525             drawOpt = isFirst ? "e" : "esames";
526             histo->SetLineColor(srcColors[isrc-1]);
527             histo->SetMarkerColor(srcColors[isrc-1]);
528             histo->SetMarkerStyle(20+4*icharge);
529             histo->Draw(drawOpt.Data());
530             TPaveStats* paveStats = (TPaveStats*)histo->FindObject("stats");
531             if ( paveStats ) paveStats->SetTextColor(srcColors[isrc-1]);
532             if ( iproj == 0 ) {
533               TString legEntry = ( itype == 0 ) ? fChargeKeys->At(icharge)->GetName() : "";
534               if ( isMC ) legEntry += Form(" %s", srcAxis->GetBinLabel(isrc));
535               if ( ! legEntry.IsNull() ) legKine->AddEntry(histo,legEntry.Data(), "lp");
536             }
537           } // loop on mu charge
538         } // loop on track sources
539       } // loop on projections
540       if ( legKine && legKine->GetNRows() > 0 ) {
541         canKine->cd(1);
542         legKine->Draw("same");
543       }
544     } // loop on grid sparse
545   } // loop on types
546   
547   
548   for ( Int_t igrid=0; igrid<3; igrid++ ) {
549     if ( gridSparseArray[igrid]->IsA() == AliCFEffGrid::Class() ) continue;
550     SetSparseRange(gridSparseArray[igrid], kHvarCharge, "", 1, gridSparseArray[igrid]->GetAxis(kHvarCharge)->GetNbins(), "USEBIN"); // Reset range
551   } // loop on container steps
552   
553   //////////////////////
554   // Event statistics //
555   //////////////////////
556   printf("\nTotal analyzed events:\n");
557   TString evtSel = Form("trigger:%s", trigClassName.Data());
558   fEventCounters->PrintSum(evtSel.Data());
559   printf("Physics selected analyzed events:\n");
560   evtSel = Form("trigger:%s/selected:yes", trigClassName.Data());
561   fEventCounters->PrintSum(evtSel.Data());
562   
563   TString countPhysSel = "any";
564   if ( physSel.Contains(fPhysSelKeys->At(kPhysSelPass)->GetName()) ) countPhysSel = "yes";
565   else if ( physSel.Contains(fPhysSelKeys->At(kPhysSelReject)->GetName()) ) countPhysSel="no";
566   countPhysSel.Prepend("selected:");
567   printf("Analyzed events vs. centrality:\n");
568   evtSel = Form("trigger:%s/%s", trigClassName.Data(), countPhysSel.Data());
569   fEventCounters->Print("centrality",evtSel.Data(),kTRUE);
570   
571   
572   TString outFilename = Form("/tmp/out%s.root", GetName());
573   TFile* outFile = new TFile(outFilename.Data(),"RECREATE");
574   for ( Int_t icharge=0; icharge<3; icharge++ ) {
575     histoList[icharge].Write();
576   }
577   histoListRatio.Write();
578   outFile->Close();
579   printf("\nCreating file %s\n", outFilename.Data(
580          ));
581   
582   ///////////////////
583   // Vertex method //
584   ///////////////////
585   if ( ! furtherOpt.Contains("VERTEX") ) return;
586   igroup1++;
587   TH1* eventVertex = (TH1*)GetSum(physSel, minBiasTrig, centralityRange, "hIpVtx");
588   if ( ! eventVertex ) return;
589   Double_t minZ = -9.99, maxZ = 9.99;
590   Double_t meanZ = 0., sigmaZ = 4.;
591   Double_t nSigma = 2.;
592   TString fitOpt = "R0S";
593   Bool_t fixFitRange = kFALSE;
594   TString fitFormula = Form("[0]+[1]*(x+[2])");
595   
596   // Get vertex shape
597   if ( eventVertex->GetSumw2N() == 0 ) eventVertex->Sumw2();
598   Double_t eventVtxIntegral = eventVertex->Integral(0,eventVertex->GetNbinsX()+1); // Include under/overflow
599   printf("Event vertex integral %.0f\n\n", eventVtxIntegral);
600   if ( eventVtxIntegral <= 0. ) return;
601   eventVertex->Scale(1./eventVtxIntegral);
602   printf("\nFit MB vertex\n");
603   eventVertex->Fit("gaus",fitOpt.Data(),"",minZ,maxZ);
604   TF1* vtxFit = (TF1*)eventVertex->GetListOfFunctions()->FindObject("gaus");
605   currName = "vtxIntegrated";
606   can = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
607   can->SetLogy();
608   eventVertex->Draw();
609   vtxFit->Draw("same");
610   
611   
612   enum {kRecoHF, kRecoBkg, kInputHF, kInputDecay, kRecoAll, kNrecoHistos};
613   TString baseRecoName[kNrecoHistos] = {"RecoHF", "RecoBkg", "InputHF", "InputDecay", "RecoAll"};
614   TArrayI sumMothers[kNrecoHistos];
615   sumMothers[kRecoHF].Set(0);
616   sumMothers[kRecoBkg].Set(0);
617   sumMothers[kInputHF].Set(3);
618   
619   sumMothers[kInputHF][0] = srcAxis->FindBin(fSrcKeys->At(kCharmMu)->GetName());
620   sumMothers[kInputHF][1] = srcAxis->FindBin(fSrcKeys->At(kBeautyMu)->GetName());
621   sumMothers[kInputHF][2] = srcAxis->FindBin(fSrcKeys->At(kQuarkoniumMu)->GetName());
622   sumMothers[kInputDecay].Set(1);
623   sumMothers[kInputDecay][0] = srcAxis->FindBin(fSrcKeys->At(kDecayMu)->GetName());
624   sumMothers[kRecoAll].Set(srcAxis->GetNbins());
625   for ( Int_t isrc=0; isrc<srcAxis->GetNbins(); ++isrc ) {
626     sumMothers[kRecoAll][isrc] = isrc+1;
627   }
628   
629   meanZ = vtxFit->GetParameter(1);
630   sigmaZ = vtxFit->GetParameter(2);
631   
632   Double_t minZfit = ( fixFitRange ) ? minZ : meanZ - nSigma*sigmaZ;
633   Double_t maxZfit = ( fixFitRange ) ? maxZ : meanZ + nSigma*sigmaZ;
634   
635   TF1* fitFunc = new TF1("fitFunc", fitFormula.Data(), minZ, maxZ);
636   fitFunc->SetLineColor(2);
637   fitFunc->SetParNames("Line norm", "Line slope", "Free path");
638   const Double_t kFreePath = 153.; // 150.; // 130.; // cm
639   //fitFunc->SetParameters(0.,1.);
640   fitFunc->FixParameter(2, kFreePath);
641   
642   AliCFGridSparse* gridSparse = cfContainer->GetGrid(kStepReconstructed);
643   TAxis* ptAxis = gridSparse->GetAxis(kHvarPt);
644   
645   Double_t slope = 0.;
646   Double_t limitNorm = 0., limitSlope = 0.;
647   Int_t firstPtBin = 0, lastPtBin = 0;
648   
649   gStyle->SetOptFit(1111);
650   
651   for ( Int_t itheta=0; itheta<kNthetaAbs; ++itheta ) {
652     igroup2++;
653     SetSparseRange(gridSparse, kHvarThetaAbs, "", itheta+1, itheta+1, "USEBIN");
654     SetSparseRange(gridSparse, kHvarPt, "", 1, ptAxis->GetNbins(), "USEBIN");
655     TH1* recoHisto[kNrecoHistos];
656     for ( Int_t ireco=0; ireco<kNrecoHistos; ++ireco ) {
657       recoHisto[ireco] = gridSparse->Project(kHvarPt);
658       histoName = Form("%sMuon_%s", baseRecoName[ireco].Data(), fThetaAbsKeys->At(itheta)->GetName());
659       recoHisto[ireco]->SetName(histoName.Data());
660       recoHisto[ireco]->SetTitle(histoName.Data());
661       recoHisto[ireco]->Reset();
662       recoHisto[ireco]->Sumw2();
663       for ( Int_t isrc=0; isrc<sumMothers[ireco].GetSize(); ++isrc ) {
664         SetSparseRange(gridSparse, kHvarMotherType, "", sumMothers[ireco][isrc], sumMothers[ireco][isrc], "USEBIN");
665         TH1* auxHisto = gridSparse->Project(kHvarPt);
666         recoHisto[ireco]->Add(auxHisto);
667         delete auxHisto;
668       }
669     }
670     SetSparseRange(gridSparse, kHvarMotherType, "", firstPtBin, lastSrcBin, "USEBIN");
671     Int_t currDraw = 0;
672     
673     for ( Int_t ibinpt=0; ibinpt<=ptAxis->GetNbins(); ++ibinpt ) {
674       firstPtBin = ibinpt;
675       lastPtBin = ( ibinpt == 0 ) ? ptAxis->GetNbins() : ibinpt;
676       SetSparseRange(gridSparse, kHvarPt, "", firstPtBin, lastPtBin, "USEBIN");
677       TH1* histo = gridSparse->Project(kHvarVz);
678       histo->SetName(Form("hVtx_%s_%s_ptBin%i", cfContainer->GetStepTitle(kStepReconstructed), fThetaAbsKeys->At(itheta)->GetName(), ibinpt));
679       if ( histo->Integral() < 100. ) break;
680       printf("\nFit %.2f < pt < %.2f (entries %g)\n", ptAxis->GetBinLowEdge(firstPtBin), ptAxis->GetBinUpEdge(lastPtBin), histo->GetEntries());
681       histo->Divide(eventVertex);
682       Double_t norm = histo->GetBinContent(histo->FindBin(0.));
683       histo->GetYaxis()->SetTitle("#frac{dN_{#mu}}{dv_{z}} / #left(#frac{1}{N_{MB}}#frac{dN_{MB}}{dv_{z}}#right)");
684       histo->SetTitle(Form("%g < p_{T} (GeV/c) < %g",ptAxis->GetBinLowEdge(firstPtBin), ptAxis->GetBinUpEdge(lastPtBin)));
685       slope = ( histo->GetBinContent(histo->FindBin(meanZ+sigmaZ)) -
686                histo->GetBinContent(histo->FindBin(meanZ-sigmaZ)) ) / ( 2. * sigmaZ );
687       
688       if ( slope < 0. ) slope = norm/kFreePath;
689       
690       // Try to fit twice: it fit fails the first time
691       // set some limits on parameters
692       for ( Int_t itry=0; itry<2; itry++ ) {
693         fitFunc->SetParameter(0, norm);
694         fitFunc->SetParameter(1, slope);
695         if ( itry > 0 ) {
696           limitNorm = 2.*histo->Integral();
697           limitSlope = 2.*histo->Integral()/kFreePath;
698           //fitFunc->SetParLimits(0, 0., limitNorm); // REMEMBER TO CHECK
699           fitFunc->SetParLimits(1, 0., limitSlope); // REMEMBER TO CHECK
700           printf("Norm 0. < %f < %f  slope  0. < %f < %f\n", norm, limitNorm, slope, limitSlope);
701         }
702         TFitResultPtr fitRes = histo->Fit(fitFunc, fitOpt.Data(), "", minZfit, maxZfit);
703         
704         //      if ( gMinuit->fCstatu.Contains("CONVERGED") &&
705         if ( ((Int_t)fitRes) == 0 &&
706             fitFunc->GetParameter(0) > 0. &&
707             fitFunc->GetParameter(1) > 0. )
708           break;
709         else if ( furtherOpt.Contains("REFIT") ) printf("Re-fit with limits\n");
710         else {
711           printf("Warning: fit problems !!!\n");
712           break;
713         }
714       }
715       
716       Double_t p0 = fitFunc->GetParameter(0);
717       Double_t p0err = fitFunc->GetParError(0);
718       Double_t p1 = fitFunc->GetParameter(1);
719       Double_t p1err = fitFunc->GetParError(1);
720       
721       Double_t nullVz = ( p1 != 0. ) ? -p0/p1 : 0.;
722       Double_t nullVzErr = ( p0 != 0. && p1 != 0. ) ? TMath::Abs(nullVz) * TMath::Sqrt(p0err*p0err/(p0*p0) + p1err*p1err/(p1*p1) ) : 0.;
723       
724       printf("Null value at %f +- %f\n", nullVz - kFreePath, nullVzErr);
725       
726       recoHisto[kRecoHF]->SetBinContent(ibinpt, p0);
727       recoHisto[kRecoHF]->SetBinError(ibinpt, p0err);
728       recoHisto[kRecoBkg]->SetBinContent(ibinpt, ( kFreePath + meanZ ) * p1);
729       recoHisto[kRecoBkg]->SetBinError(ibinpt, ( kFreePath + meanZ ) * p1err);
730       if ( currDraw%4 == 0 ){
731         currName = Form("vtx_%s_PtBin%i",fThetaAbsKeys->At(itheta)->GetName(), ibinpt);
732         can = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
733         can->Divide(2,2);
734       }
735       can->cd( currDraw%4 + 1 );
736       can->SetLogy();
737       histo->Draw();
738       fitFunc->DrawCopy("same");
739       currDraw++;
740     } // loop on pt bins
741     SetSparseRange(gridSparse, kHvarMotherType, "", firstSrcBin, lastSrcBin, "USEBIN");
742     currName = Form("recoPt_%s",fThetaAbsKeys->At(itheta)->GetName());
743     can = new TCanvas(currName.Data(),currName.Data(),(igroup1+1)*xshift,igroup2*yshift,600,600);
744     TLegend* leg = new TLegend(0.6, 0.6, 0.8, 0.8);
745     drawOpt = "e";
746     for ( Int_t ireco=0; ireco<kNrecoHistos-1; ++ireco ) {
747       if ( recoHisto[ireco]->GetEntries() == 0. ) continue;
748       TH1* ratio = (TH1*)recoHisto[ireco]->Clone(Form("%s_ratio", recoHisto[ireco]->GetName()));
749       ratio->Divide(recoHisto[kRecoAll]);
750       ratio->SetLineColor(srcColors[ireco]);
751       ratio->SetMarkerColor(srcColors[ireco]);
752       ratio->SetMarkerStyle(20+ireco);
753       ratio->GetYaxis()->SetTitle("fraction of total");
754       ratio->Draw(drawOpt.Data());
755       leg->AddEntry(ratio,baseRecoName[ireco].Data(), "lp");
756       drawOpt = "esame";
757     }
758     leg->Draw("same");
759   } // loop on theta abs
760 }