1 /**************************************************************************
2 * Copyright(c) 1998-2007, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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.
25 /// \author Diego Stocco
26 //-----------------------------------------------------------------------------
28 #define AliAnalysisTaskSingleMu_cxx
30 #include "AliAnalysisTaskSingleMu.h"
40 #include "TObjString.h"
41 #include "TObjArray.h"
44 //#include "TMCProcess.h"
47 #include "TPaveStats.h"
48 #include "TFitResultPtr.h"
50 #include "THashList.h"
53 #include "AliAODEvent.h"
54 #include "AliAODTrack.h"
55 #include "AliAODMCParticle.h"
56 #include "AliMCEvent.h"
57 #include "AliMCParticle.h"
58 #include "AliESDEvent.h"
59 #include "AliESDMuonTrack.h"
60 #include "AliVHeader.h"
61 #include "AliAODMCHeader.h"
65 #include "AliAnalysisManager.h"
68 #include "AliCFContainer.h"
69 #include "AliCFGridSparse.h"
70 #include "AliCFEffGrid.h"
73 #include "AliVAnalysisMuon.h"
74 #include "AliMergeableCollection.h"
75 #include "AliCounterCollection.h"
76 #include "AliMuonEventCuts.h"
77 #include "AliMuonTrackCuts.h"
78 #include "AliAnalysisMuonUtility.h"
82 ClassImp(AliAnalysisTaskSingleMu) // Class implementation in ROOT context
86 //________________________________________________________________________
87 AliAnalysisTaskSingleMu::AliAnalysisTaskSingleMu() :
95 //________________________________________________________________________
96 AliAnalysisTaskSingleMu::AliAnalysisTaskSingleMu(const char *name, const AliMuonTrackCuts& cuts) :
97 AliVAnalysisMuon(name, cuts),
104 TString thetaAbsKeys = "ThetaAbs23 ThetaAbs310";
105 fThetaAbsKeys = thetaAbsKeys.Tokenize(" ");
109 //________________________________________________________________________
110 AliAnalysisTaskSingleMu::~AliAnalysisTaskSingleMu()
116 delete fThetaAbsKeys;
119 //___________________________________________________________________________
120 void AliAnalysisTaskSingleMu::MyUserCreateOutputObjects()
124 TString histoName = "", histoTitle = "";
127 Double_t vzMin = -20., vzMax = 20.;
128 TString vzName("Vz"), vzTitle("Vz"), vzUnits("cm");
130 histoName = "hIpVtx";
131 histo = new TH1F(histoName.Data(), histoName.Data(), nVzBins, vzMin, vzMax);
132 histo->SetXTitle("v_{z} (cm)");
133 AddObjectToCollection(histo, kIPVz);
137 Double_t ptMin = 0., ptMax = 80.;
138 TString ptName("Pt"), ptTitle("p_{T}"), ptUnits("GeV/c");
141 Double_t etaMin = -4.5, etaMax = -2.;
142 TString etaName("Eta"), etaTitle("#eta"), etaUnits("");
145 Double_t phiMin = 0.; Double_t phiMax = 2.*TMath::Pi();
146 TString phiName("Phi"), phiTitle("#phi"), phiUnits("rad");
148 Int_t nChargeBins = 2;
149 Double_t chargeMin = -2, chargeMax = 2.;
150 TString chargeName("Charge"), chargeTitle("charge"), chargeUnits("e");
152 Int_t nThetaAbsEndBins = 2;
153 Double_t thetaAbsEndMin = -0.5, thetaAbsEndMax = 1.5;
154 TString thetaAbsEndName("ThetaAbsEnd"), thetaAbsEndTitle("#theta_{abs}"), thetaAbsEndUnits("a.u.");
156 Int_t nMotherTypeBins = kNtrackSources;
157 Double_t motherTypeMin = -0.5, motherTypeMax = (Double_t)kNtrackSources - 0.5;
158 TString motherType("MotherType"), motherTypeTitle("motherType"), motherTypeUnits("");
160 Int_t nbins[kNvars] = {nPtBins, nEtaBins, nPhiBins, nVzBins, nChargeBins, nThetaAbsEndBins, nMotherTypeBins};
161 Double_t xmin[kNvars] = {ptMin, etaMin, phiMin, vzMin, chargeMin, thetaAbsEndMin, motherTypeMin};
162 Double_t xmax[kNvars] = {ptMax, etaMax, phiMax, vzMax, chargeMax, thetaAbsEndMax, motherTypeMax};
163 TString axisTitle[kNvars] = {ptTitle, etaTitle, phiTitle, vzTitle, chargeTitle, thetaAbsEndTitle, motherTypeTitle};
164 TString axisUnits[kNvars] = {ptUnits, etaUnits, phiUnits, vzUnits, chargeUnits, thetaAbsEndUnits, motherTypeUnits};
166 AliCFContainer* cfContainer = new AliCFContainer("SingleMuContainer","Container for tracks",kNsteps,kNvars,nbins);
168 for ( Int_t idim = 0; idim<kNvars; idim++){
169 histoTitle = Form("%s (%s)", axisTitle[idim].Data(), axisUnits[idim].Data());
170 histoTitle.ReplaceAll("()","");
172 cfContainer->SetVarTitle(idim, histoTitle.Data());
173 cfContainer->SetBinLimits(idim, xmin[idim], xmax[idim]);
176 TString stepTitle[kNsteps] = {"reconstructed", "generated"};
178 TAxis* currAxis = 0x0;
179 for (Int_t istep=0; istep<kNsteps; istep++){
180 cfContainer->SetStepTitle(istep, stepTitle[istep].Data());
181 AliCFGridSparse* gridSparse = cfContainer->GetGrid(istep);
183 currAxis = gridSparse->GetAxis(kHvarMotherType);
184 for ( Int_t ibin=0; ibin<fSrcKeys->GetEntries(); ibin++ ) {
185 currAxis->SetBinLabel(ibin+1, fSrcKeys->At(ibin)->GetName());
189 AddObjectToCollection(cfContainer, kTrackContainer);
191 histoName = "hRecoDimu";
192 TH2* histoDimu = new TH2F(histoName.Data(), histoName.Data(), 24, 0., 120., 30, 0., 150.);
193 histoDimu->SetXTitle("min. p_{T}^{#mu} (GeV/c)");
194 histoDimu->SetYTitle("M_{#mu#mu} (GeV/c^{2})");
195 AddObjectToCollection(histoDimu, kNobjectTypes);
197 histoName = "hGeneratedZ";
198 TH2* histoGenZ = static_cast<TH2*>(histoDimu->Clone(histoName.Data()));
199 histoGenZ->SetTitle(histoName.Data());
200 AddObjectToCollection(histoGenZ, kNobjectTypes+1);
203 histoName = "hZmuEtaCorr";
204 histoDimu = new TH2F(histoName.Data(), histoName.Data(), 160, -8., 8., 160, -8., 8.);
205 histoDimu->SetXTitle("#eta_{#mu-}");
206 histoDimu->SetYTitle("#eta_{#mu+}");
207 AddObjectToCollection(histoDimu, kNobjectTypes+2);
209 fMuonTrackCuts->Print("mask");
211 AliInfo(Form("Apply cut on dimuon (60<M_mumu<120 GeV/c^2) to reject Z contribution: %i", fCutOnDimu));
214 //________________________________________________________________________
215 void AliAnalysisTaskSingleMu::ProcessEvent(TString physSel, const TObjArray& selectTrigClasses, TString centrality)
218 /// Fill output objects
221 Double_t ipVz = AliAnalysisMuonUtility::GetVertexSPD(InputEvent())->GetZ();
222 Double_t ipVzMC = MCEvent() ? AliAnalysisMuonUtility::GetMCVertexZ(InputEvent(),MCEvent()) : 0.;
224 for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
225 TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
226 ((TH1*)GetMergeableObject(physSel, trigClassName, centrality, "hIpVtx"))->Fill(ipVz);
229 // 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
230 // if ( isPileupFromSPD ) return; // UNCOMMENT TO REJECT PILEUP
232 Double_t containerInput[kNvars];
233 AliVParticle* track = 0x0;
235 Int_t nSteps = MCEvent() ? 2 : 1;
236 for ( Int_t istep = 0; istep<nSteps; ++istep ) {
237 Int_t nTracks = ( istep == kStepReconstructed ) ? AliAnalysisMuonUtility::GetNTracks(InputEvent()) : MCEvent()->GetNumberOfTracks();
239 TObjArray selectedTracks(nTracks);
240 TArrayI trackSources(nTracks);
241 TArrayD trackWgt(nTracks);
243 for (Int_t itrack = 0; itrack < nTracks; itrack++) {
244 track = ( istep == kStepReconstructed ) ? AliAnalysisMuonUtility::GetTrack(itrack,InputEvent()) : MCEvent()->GetTrack(itrack);
246 // In case of MC we usually ask that the particle is a muon
247 // However, in W or Z simulations, Pythia stores both the initial muon
248 // (before ISR, FSR and kt kick) and the final state one.
249 // The first muon is of course there only for information and should be rejected.
250 // The Pythia code for initial state particles is 21
251 Bool_t isSelected = ( istep == kStepReconstructed ) ? fMuonTrackCuts->IsSelected(track) : ( TMath::Abs(track->PdgCode()) == 13 && AliAnalysisMuonUtility::GetStatusCode(track) != 21 );
252 if ( ! isSelected ) continue;
254 selectedTracks.AddAt(track,itrack);
255 trackSources[itrack] = GetParticleType(track);
257 TObject* wgtObj = GetWeight(fSrcKeys->At(trackSources[itrack])->GetName());
260 AliVParticle* mcTrack = ( istep == kStepReconstructed ) ? MCEvent()->GetTrack(track->GetLabel()) : track;
261 if ( wgtObj->IsA() == TF1::Class() ) trackWgt[itrack] = static_cast<TF1*>(wgtObj)->Eval(mcTrack->Pt());
262 else if ( wgtObj->IsA()->InheritsFrom(TH1::Class()) ) {
263 TH1* wgtHisto = static_cast<TH1*>(wgtObj);
264 trackWgt[itrack] = wgtHisto->GetBinContent(wgtHisto->GetXaxis()->FindBin(mcTrack->Pt()));
266 AliDebug(3,Form("Apply weights %s: pt %g gen pt %g weight %g",wgtObj->GetName(),track->Pt(),mcTrack->Pt(),trackWgt[itrack]));
267 // Int_t iancestor = fUtilityMuonAncestor->GetAncestor(track,MCEvent());
268 // AliVParticle* motherPart = MCEvent()->GetTrack(iancestor);
269 // trackWgt[itrack] = beautyMuWgt->GetBinContent(beautyMuWgt->GetXaxis()->FindBin(motherPart->Pt()));
273 // Loop on selected tracks
274 TArrayI rejectTrack(nTracks);
275 for ( Int_t itrack=0; itrack<nTracks; itrack++) {
276 track = static_cast<AliVParticle*>(selectedTracks.At(itrack));
277 if ( ! track ) continue;
280 for ( Int_t jtrack=itrack+1; jtrack<nTracks; jtrack++ ) {
281 AliVParticle* auxTrack = static_cast<AliVParticle*>(selectedTracks.At(jtrack));
282 if ( ! auxTrack ) continue;
283 if ( track->Charge() * auxTrack->Charge() >= 0 ) continue;
285 TLorentzVector dimuPair = AliAnalysisMuonUtility::GetTrackPair(track,auxTrack);
286 Double_t ptMin = TMath::Min(track->Pt(),auxTrack->Pt());
287 Double_t invMass = dimuPair.M();
288 if ( fCutOnDimu && invMass > 60. && invMass < 120. ) {
289 rejectTrack[itrack] = 1;
290 rejectTrack[jtrack] = 1;
292 Double_t dimuWgt = trackWgt[itrack] * trackWgt[jtrack];
293 if ( istep == kStepReconstructed ) {
294 for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
295 TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
296 if ( ! fMuonTrackCuts->TrackPtCutMatchTrigClass(track, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) || ! fMuonTrackCuts->TrackPtCutMatchTrigClass(auxTrack, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) ) continue;
297 ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hRecoDimu"))->Fill(ptMin,invMass,dimuWgt);
298 } // loop on triggers
301 if ( trackSources[itrack] == kZbosonMu && trackSources[jtrack] == kZbosonMu ) {
302 Bool_t isAccepted = kTRUE;
303 AliVParticle* muDaughter[2] = {0x0, 0x0};
304 for ( Int_t imu=0; imu<2; imu++ ) {
305 AliVParticle* currPart = ( imu == 0 ) ? track : auxTrack;
306 if ( currPart->Charge() < 0. ) muDaughter[0] = currPart;
307 else muDaughter[1] = currPart;
308 if ( currPart->Eta() < -4.5 || currPart->Eta() > -2. ) {
311 } // loop on muons in the pair
313 Double_t pairRapidity = dimuPair.Rapidity();
314 if ( pairRapidity < -4. || pairRapidity > -2.5 ) isAccepted = kFALSE;
315 //printf("Rapidity Z %g pair %g\n",track->Y(), pairRapidity);
317 for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
318 TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
319 if ( isAccepted ) ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hGeneratedZ"))->Fill(ptMin,invMass,dimuWgt);
320 ((TH2*)GetMergeableObject(physSel, trigClassName, centrality, "hZmuEtaCorr"))->Fill(muDaughter[0]->Eta(),muDaughter[1]->Eta(),dimuWgt);
321 } // loop on selected trig
322 } // both muons from Z
323 } // kStepGeneratedMC
324 } // loop on auxiliary tracks
325 if ( rejectTrack[itrack] > 0 ) continue;
327 Double_t thetaAbsEndDeg = 0;
328 if ( istep == kStepReconstructed ) {
329 thetaAbsEndDeg = AliAnalysisMuonUtility::GetThetaAbsDeg(track);
332 thetaAbsEndDeg = ( TMath::Pi()-track->Theta() ) * TMath::RadToDeg();
334 Int_t thetaAbsBin = ( thetaAbsEndDeg < 3. ) ? kThetaAbs23 : kThetaAbs310;
336 containerInput[kHvarPt] = track->Pt();
337 containerInput[kHvarEta] = track->Eta();
338 containerInput[kHvarPhi] = track->Phi();
339 containerInput[kHvarVz] = ( istep == kStepReconstructed ) ? ipVz : ipVzMC;
340 containerInput[kHvarCharge] = track->Charge()/3.;
341 containerInput[kHvarThetaAbs] = (Double_t)thetaAbsBin;
342 containerInput[kHvarMotherType] = (Double_t)trackSources[itrack];
344 for ( Int_t itrig=0; itrig<selectTrigClasses.GetEntries(); ++itrig ) {
345 TString trigClassName = ((TObjString*)selectTrigClasses.At(itrig))->GetString();
346 if ( istep == kStepReconstructed && ! fMuonTrackCuts->TrackPtCutMatchTrigClass(track, fMuonEventCuts->GetTrigClassPtCutLevel(trigClassName)) ) continue;
347 ((AliCFContainer*)GetMergeableObject(physSel, trigClassName, centrality, "SingleMuContainer"))->Fill(containerInput,istep,trackWgt[itrack]);
348 } // loop on selected trigger classes
350 } // loop on container steps
354 //________________________________________________________________________
355 void AliAnalysisTaskSingleMu::Terminate(Option_t *) {
357 /// Draw some histograms at the end.
360 AliVAnalysisMuon::Terminate("");
362 if ( ! fMergeableCollection ) return;
364 TString physSel = fTerminateOptions->At(0)->GetName();
365 TString trigClassName = fTerminateOptions->At(1)->GetName();
366 TString centralityRange = fTerminateOptions->At(2)->GetName();
367 TString furtherOpt = fTerminateOptions->At(3)->GetName();
369 TString minBiasTrig = "";
370 TObjArray* optArr = furtherOpt.Tokenize(" ");
371 TString currName = "";
372 for ( Int_t iopt=0; iopt<optArr->GetEntries(); iopt++ ) {
373 currName = optArr->At(iopt)->GetName();
374 if ( currName.Contains("-B-") || currName.Contains("ANY") ) minBiasTrig = currName;
378 furtherOpt.ToUpper();
379 Bool_t plotChargeAsymmetry = furtherOpt.Contains("ASYM");
381 AliCFContainer* inContainer = static_cast<AliCFContainer*> ( GetSum(physSel,trigClassName,centralityRange,"SingleMuContainer") );
382 if ( ! inContainer ) return;
384 AliCFContainer* cfContainer = inContainer;
386 if ( ! furtherOpt.Contains("GENINTRIGCLASS") && trigClassName != "ANY" ) {
387 // The output container contains both the reconstructed and (in MC)
388 // the generated muons in a specific trigger class.
389 // Since the trigger pt level of the track is matched to the trigger class,
390 // analyzing the MUHigh trigger (for example) is equivalent of analyzing
392 // However, in this way, the generated muons distribution depend
393 // on a condition on the reconstructed track.
394 // If we calculate the Acc.xEff. as the ratio of reconstructed and
395 // generated muons IN A TRIGGER CLASS, we are biasing the final result.
396 // The correct value of Acc.xEff. is hence obtained as the distribution
397 // of reconstructed tracks in a muon trigger class, divided by the
398 // total number of generated class (which is in the "class" ANY).
399 // The following code sets the generated muons as the one in the class ANY.
400 // The feature is the default. If you want the generated muons in the same
401 // trigger class as the generated tracks, use option "GENINTRIGCLASS"
402 AliCFContainer* fullMCcontainer = static_cast<AliCFContainer*> ( GetSum(physSel,"ANY",centralityRange,"SingleMuContainer") );
403 if ( fullMCcontainer ) {
404 cfContainer = static_cast<AliCFContainer*>(cfContainer->Clone("SingleMuContainerCombo"));
405 cfContainer->SetGrid(kStepGeneratedMC,fullMCcontainer->GetGrid(kStepGeneratedMC));
409 AliCFEffGrid* effSparse = new AliCFEffGrid(Form("eff%s", cfContainer->GetName()),Form("Efficiency %s", cfContainer->GetTitle()),*cfContainer);
410 effSparse->CalculateEfficiency(kStepReconstructed, kStepGeneratedMC);
412 AliCFGridSparse* gridSparseArray[3] = {effSparse->GetNum(), effSparse->GetDen(), effSparse};
413 TString gridSparseName[3] = {cfContainer->GetStepTitle(kStepReconstructed), cfContainer->GetStepTitle(kStepGeneratedMC), "Efficiency"};
415 Int_t srcColors[kNtrackSources] = {kBlack, kRed, kSpring, kTeal, kBlue, kViolet, kMagenta, kOrange, kGray};
416 // TString allSrcNames = "";
417 // for ( Int_t isrc=0; isrc<kNtrackSources; ++isrc ) {
418 // if ( ! allSrcNames.IsNull() ) allSrcNames.Append(" ");
419 // allSrcNames += fSrcKeys->At(isrc)->GetName();
428 Bool_t isMC = furtherOpt.Contains("MC");
430 TAxis* srcAxis = gridSparseArray[0]->GetAxis(kHvarMotherType);
431 Int_t unIdBin = srcAxis->FindBin(fSrcKeys->At(kUnidentified)->GetName());
432 if ( unIdBin < 1 ) unIdBin = srcAxis->GetNbins();
434 Int_t firstSrcBin = ( isMC ) ? 1 : unIdBin;
435 Int_t lastSrcBin = ( isMC ) ? srcAxis->GetNbins() - 1 : unIdBin;
436 if ( ! isMC ) srcColors[unIdBin-1] = 1;
441 TString chargeNames[3] = {fChargeKeys->At(0)->GetName(), fChargeKeys->At(1)->GetName(), "Total"};
442 THashList histoList[3];
443 for ( Int_t icharge=0; icharge<3; icharge++ ) {
444 histoList[icharge].SetName(chargeNames[icharge].Data());
446 for ( Int_t isrc = firstSrcBin; isrc <= lastSrcBin; ++isrc ) {
447 for ( Int_t icharge=0; icharge<3; ++icharge ) {
448 Int_t icharge1 = ( icharge < 2 ) ? icharge : 0;
449 Int_t icharge2 = ( icharge < 2 ) ? icharge : 1;
450 for ( Int_t igrid=0; igrid<3; ++igrid ) {
451 if ( gridSparseArray[igrid]->GetEntries() == 0. ) break;
452 if ( gridSparseArray[igrid]->IsA() != AliCFEffGrid::Class() ) {
453 SetSparseRange(gridSparseArray[igrid], kHvarEta, "", -3.999, -2.501);
454 SetSparseRange(gridSparseArray[igrid], kHvarMotherType, "", isrc, isrc, "USEBIN");
455 SetSparseRange(gridSparseArray[igrid], kHvarCharge, "", icharge1+1, icharge2+1, "USEBIN");
457 TH1* histo = gridSparseArray[igrid]->Project(kHvarPt, kHvarEta);
458 histo->SetName(Form("hPtEta_%s_%s_%s", gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), chargeNames[icharge].Data()));
459 if ( histo->Integral() > 0 ) histoList[icharge].Add(histo);
460 for ( Int_t iproj=0; iproj<4; ++iproj ) {
461 histo = gridSparseArray[igrid]->Project(iproj);
462 histo->SetName(Form("proj%i_%s_%s_%s", iproj, gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), chargeNames[icharge].Data()));
463 if ( histo->Integral() > 0 ) histoList[icharge].Add(histo);
464 } // loop on projections
465 } // loop on grid sparse
467 } // loop on track sources
469 // Get charge asymmetry or mu+/mu-
470 THashList histoListRatio;
471 TString basePlotName = plotChargeAsymmetry ? "ChargeAsymmetry" : "ChargeRatio";
472 histoListRatio.SetName(basePlotName.Data());
473 Int_t baseCharge = 1;
474 Int_t auxCharge = 1-baseCharge;
475 for ( Int_t ihisto=0; ihisto<histoList[baseCharge].GetEntries(); ihisto++ ) {
476 TObject* obj = histoList[baseCharge].At(ihisto);
477 TString histoName = obj->GetName();
478 if ( histoName.Contains(gridSparseName[2].Data()) ) continue;
479 TString searchName = histoName;
480 searchName.ReplaceAll(fChargeKeys->At(baseCharge)->GetName(), fChargeKeys->At(auxCharge)->GetName());
481 TH1* auxHisto = static_cast<TH1*> (histoList[auxCharge].FindObject(searchName.Data()));
482 if ( ! auxHisto ) continue;
483 histoName.ReplaceAll(fChargeKeys->At(baseCharge)->GetName(),basePlotName.Data());
484 TH1* histo = static_cast<TH1*> (obj->Clone(histoName.Data()));
485 if ( plotChargeAsymmetry ) {
486 histo->Add(auxHisto, -1.);
487 // h2 + h1 = 2xh2 + (h1-h2)
488 auxHisto->Add(auxHisto, histo, 2.);
490 histo->Divide(auxHisto);
491 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());
492 axisTitle.ReplaceAll("MuPlus","#mu^{+}");
493 axisTitle.ReplaceAll("MuMinus","#mu^{-}");
494 histo->GetYaxis()->SetTitle(axisTitle.Data());
495 histo->SetStats(kFALSE);
496 histoListRatio.Add(histo);
500 TString histoName = "", drawOpt = "";
501 for ( Int_t itype=0; itype<3; itype++ ) {
502 THashList* currList = 0x0;
504 if ( itype == 1 ) currList = &histoListRatio;
505 else if ( itype == 2 ) currList = &histoList[2];
507 for ( Int_t igrid=0; igrid<3; ++igrid ) {
509 TCanvas* canKine = 0x0;
510 TLegend* legKine = 0x0;
511 for ( Int_t iproj=0; iproj<4; ++iproj ) {
512 for ( Int_t isrc = firstSrcBin; isrc <= lastSrcBin; ++isrc ) {
513 for ( Int_t icharge=0; icharge<nCharges; ++icharge ) {
514 if ( itype == 0 ) currList = &histoList[icharge];
515 histoName = Form("proj%i_%s_%s_%s", iproj, gridSparseName[igrid].Data(), srcAxis->GetBinLabel(isrc), currList->GetName());
516 TH1* histo = static_cast<TH1*>(currList->FindObject(histoName.Data()));
517 if ( ! histo ) continue;
521 currName = Form("%s_%s_%s", GetName(), currList->GetName(), gridSparseName[igrid].Data());
522 canKine = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
523 canKine->Divide(2,2);
524 legKine = new TLegend(0.6, 0.6, 0.8, 0.8);
526 canKine->cd(iproj+1);
528 if ( ( iproj == kHvarPt || iproj == kHvarVz ) && gridSparseArray[igrid]->IsA() != AliCFEffGrid::Class() ) gPad->SetLogy();
530 Bool_t isFirst = ( gPad->GetListOfPrimitives()->GetEntries() == 0 );
531 drawOpt = isFirst ? "e" : "esames";
532 histo->SetLineColor(srcColors[isrc-1]);
533 histo->SetMarkerColor(srcColors[isrc-1]);
534 histo->SetMarkerStyle(20+4*icharge);
535 histo->Draw(drawOpt.Data());
536 TPaveStats* paveStats = (TPaveStats*)histo->FindObject("stats");
537 if ( paveStats ) paveStats->SetTextColor(srcColors[isrc-1]);
539 TString legEntry = ( itype == 0 ) ? fChargeKeys->At(icharge)->GetName() : "";
540 if ( isMC ) legEntry += Form(" %s", srcAxis->GetBinLabel(isrc));
541 if ( ! legEntry.IsNull() ) legKine->AddEntry(histo,legEntry.Data(), "lp");
543 } // loop on mu charge
544 } // loop on track sources
545 } // loop on projections
546 if ( legKine && legKine->GetNRows() > 0 ) {
548 legKine->Draw("same");
550 } // loop on grid sparse
554 for ( Int_t igrid=0; igrid<3; igrid++ ) {
555 if ( gridSparseArray[igrid]->IsA() == AliCFEffGrid::Class() ) continue;
556 SetSparseRange(gridSparseArray[igrid], kHvarCharge, "", 1, gridSparseArray[igrid]->GetAxis(kHvarCharge)->GetNbins(), "USEBIN"); // Reset range
557 } // loop on container steps
559 //////////////////////
560 // Event statistics //
561 //////////////////////
562 printf("\nTotal analyzed events:\n");
563 TString evtSel = Form("trigger:%s", trigClassName.Data());
564 fEventCounters->PrintSum(evtSel.Data());
565 printf("Physics selected analyzed events:\n");
566 evtSel = Form("trigger:%s/selected:yes", trigClassName.Data());
567 fEventCounters->PrintSum(evtSel.Data());
569 TString countPhysSel = "any";
570 if ( physSel.Contains(fPhysSelKeys->At(kPhysSelPass)->GetName()) ) countPhysSel = "yes";
571 else if ( physSel.Contains(fPhysSelKeys->At(kPhysSelReject)->GetName()) ) countPhysSel="no";
572 countPhysSel.Prepend("selected:");
573 printf("Analyzed events vs. centrality:\n");
574 evtSel = Form("trigger:%s/%s", trigClassName.Data(), countPhysSel.Data());
575 fEventCounters->Print("centrality",evtSel.Data(),kTRUE);
578 TString outFilename = Form("/tmp/out%s.root", GetName());
579 TFile* outFile = new TFile(outFilename.Data(),"RECREATE");
580 for ( Int_t icharge=0; icharge<3; icharge++ ) {
581 histoList[icharge].Write();
583 histoListRatio.Write();
585 printf("\nCreating file %s\n", outFilename.Data(
591 if ( ! furtherOpt.Contains("VERTEX") ) return;
593 TH1* eventVertex = (TH1*)GetSum(physSel, minBiasTrig, centralityRange, "hIpVtx");
594 if ( ! eventVertex ) return;
595 Double_t minZ = -9.99, maxZ = 9.99;
596 Double_t meanZ = 0., sigmaZ = 4.;
597 Double_t nSigma = 2.;
598 TString fitOpt = "R0S";
599 Bool_t fixFitRange = kFALSE;
600 TString fitFormula = Form("[0]+[1]*(x+[2])");
603 if ( eventVertex->GetSumw2N() == 0 ) eventVertex->Sumw2();
604 Double_t eventVtxIntegral = eventVertex->Integral(0,eventVertex->GetNbinsX()+1); // Include under/overflow
605 printf("Event vertex integral %.0f\n\n", eventVtxIntegral);
606 if ( eventVtxIntegral <= 0. ) return;
607 eventVertex->Scale(1./eventVtxIntegral);
608 printf("\nFit MB vertex\n");
609 eventVertex->Fit("gaus",fitOpt.Data(),"",minZ,maxZ);
610 TF1* vtxFit = (TF1*)eventVertex->GetListOfFunctions()->FindObject("gaus");
611 currName = "vtxIntegrated";
612 can = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
615 vtxFit->Draw("same");
618 enum {kRecoHF, kRecoBkg, kInputHF, kInputDecay, kRecoAll, kNrecoHistos};
619 TString baseRecoName[kNrecoHistos] = {"RecoHF", "RecoBkg", "InputHF", "InputDecay", "RecoAll"};
620 TArrayI sumMothers[kNrecoHistos];
621 sumMothers[kRecoHF].Set(0);
622 sumMothers[kRecoBkg].Set(0);
623 sumMothers[kInputHF].Set(3);
625 sumMothers[kInputHF][0] = srcAxis->FindBin(fSrcKeys->At(kCharmMu)->GetName());
626 sumMothers[kInputHF][1] = srcAxis->FindBin(fSrcKeys->At(kBeautyMu)->GetName());
627 sumMothers[kInputHF][2] = srcAxis->FindBin(fSrcKeys->At(kQuarkoniumMu)->GetName());
628 sumMothers[kInputDecay].Set(1);
629 sumMothers[kInputDecay][0] = srcAxis->FindBin(fSrcKeys->At(kDecayMu)->GetName());
630 sumMothers[kRecoAll].Set(srcAxis->GetNbins());
631 for ( Int_t isrc=0; isrc<srcAxis->GetNbins(); ++isrc ) {
632 sumMothers[kRecoAll][isrc] = isrc+1;
635 meanZ = vtxFit->GetParameter(1);
636 sigmaZ = vtxFit->GetParameter(2);
638 Double_t minZfit = ( fixFitRange ) ? minZ : meanZ - nSigma*sigmaZ;
639 Double_t maxZfit = ( fixFitRange ) ? maxZ : meanZ + nSigma*sigmaZ;
641 TF1* fitFunc = new TF1("fitFunc", fitFormula.Data(), minZ, maxZ);
642 fitFunc->SetLineColor(2);
643 fitFunc->SetParNames("Line norm", "Line slope", "Free path");
644 const Double_t kFreePath = 153.; // 150.; // 130.; // cm
645 //fitFunc->SetParameters(0.,1.);
646 fitFunc->FixParameter(2, kFreePath);
648 AliCFGridSparse* gridSparse = cfContainer->GetGrid(kStepReconstructed);
649 TAxis* ptAxis = gridSparse->GetAxis(kHvarPt);
652 Double_t limitNorm = 0., limitSlope = 0.;
653 Int_t firstPtBin = 0, lastPtBin = 0;
655 gStyle->SetOptFit(1111);
657 for ( Int_t itheta=0; itheta<kNthetaAbs; ++itheta ) {
659 SetSparseRange(gridSparse, kHvarThetaAbs, "", itheta+1, itheta+1, "USEBIN");
660 SetSparseRange(gridSparse, kHvarPt, "", 1, ptAxis->GetNbins(), "USEBIN");
661 TH1* recoHisto[kNrecoHistos];
662 for ( Int_t ireco=0; ireco<kNrecoHistos; ++ireco ) {
663 recoHisto[ireco] = gridSparse->Project(kHvarPt);
664 histoName = Form("%sMuon_%s", baseRecoName[ireco].Data(), fThetaAbsKeys->At(itheta)->GetName());
665 recoHisto[ireco]->SetName(histoName.Data());
666 recoHisto[ireco]->SetTitle(histoName.Data());
667 recoHisto[ireco]->Reset();
668 recoHisto[ireco]->Sumw2();
669 for ( Int_t isrc=0; isrc<sumMothers[ireco].GetSize(); ++isrc ) {
670 SetSparseRange(gridSparse, kHvarMotherType, "", sumMothers[ireco][isrc], sumMothers[ireco][isrc], "USEBIN");
671 TH1* auxHisto = gridSparse->Project(kHvarPt);
672 recoHisto[ireco]->Add(auxHisto);
676 SetSparseRange(gridSparse, kHvarMotherType, "", firstPtBin, lastSrcBin, "USEBIN");
679 for ( Int_t ibinpt=0; ibinpt<=ptAxis->GetNbins(); ++ibinpt ) {
681 lastPtBin = ( ibinpt == 0 ) ? ptAxis->GetNbins() : ibinpt;
682 SetSparseRange(gridSparse, kHvarPt, "", firstPtBin, lastPtBin, "USEBIN");
683 TH1* histo = gridSparse->Project(kHvarVz);
684 histo->SetName(Form("hVtx_%s_%s_ptBin%i", cfContainer->GetStepTitle(kStepReconstructed), fThetaAbsKeys->At(itheta)->GetName(), ibinpt));
685 if ( histo->Integral() < 100. ) break;
686 printf("\nFit %.2f < pt < %.2f (entries %g)\n", ptAxis->GetBinLowEdge(firstPtBin), ptAxis->GetBinUpEdge(lastPtBin), histo->GetEntries());
687 histo->Divide(eventVertex);
688 Double_t norm = histo->GetBinContent(histo->FindBin(0.));
689 histo->GetYaxis()->SetTitle("#frac{dN_{#mu}}{dv_{z}} / #left(#frac{1}{N_{MB}}#frac{dN_{MB}}{dv_{z}}#right)");
690 histo->SetTitle(Form("%g < p_{T} (GeV/c) < %g",ptAxis->GetBinLowEdge(firstPtBin), ptAxis->GetBinUpEdge(lastPtBin)));
691 slope = ( histo->GetBinContent(histo->FindBin(meanZ+sigmaZ)) -
692 histo->GetBinContent(histo->FindBin(meanZ-sigmaZ)) ) / ( 2. * sigmaZ );
694 if ( slope < 0. ) slope = norm/kFreePath;
696 // Try to fit twice: it fit fails the first time
697 // set some limits on parameters
698 for ( Int_t itry=0; itry<2; itry++ ) {
699 fitFunc->SetParameter(0, norm);
700 fitFunc->SetParameter(1, slope);
702 limitNorm = 2.*histo->Integral();
703 limitSlope = 2.*histo->Integral()/kFreePath;
704 //fitFunc->SetParLimits(0, 0., limitNorm); // REMEMBER TO CHECK
705 fitFunc->SetParLimits(1, 0., limitSlope); // REMEMBER TO CHECK
706 printf("Norm 0. < %f < %f slope 0. < %f < %f\n", norm, limitNorm, slope, limitSlope);
708 TFitResultPtr fitRes = histo->Fit(fitFunc, fitOpt.Data(), "", minZfit, maxZfit);
710 // if ( gMinuit->fCstatu.Contains("CONVERGED") &&
711 if ( ((Int_t)fitRes) == 0 &&
712 fitFunc->GetParameter(0) > 0. &&
713 fitFunc->GetParameter(1) > 0. )
715 else if ( furtherOpt.Contains("REFIT") ) printf("Re-fit with limits\n");
717 printf("Warning: fit problems !!!\n");
722 Double_t p0 = fitFunc->GetParameter(0);
723 Double_t p0err = fitFunc->GetParError(0);
724 Double_t p1 = fitFunc->GetParameter(1);
725 Double_t p1err = fitFunc->GetParError(1);
727 Double_t nullVz = ( p1 != 0. ) ? -p0/p1 : 0.;
728 Double_t nullVzErr = ( p0 != 0. && p1 != 0. ) ? TMath::Abs(nullVz) * TMath::Sqrt(p0err*p0err/(p0*p0) + p1err*p1err/(p1*p1) ) : 0.;
730 printf("Null value at %f +- %f\n", nullVz - kFreePath, nullVzErr);
732 recoHisto[kRecoHF]->SetBinContent(ibinpt, p0);
733 recoHisto[kRecoHF]->SetBinError(ibinpt, p0err);
734 recoHisto[kRecoBkg]->SetBinContent(ibinpt, ( kFreePath + meanZ ) * p1);
735 recoHisto[kRecoBkg]->SetBinError(ibinpt, ( kFreePath + meanZ ) * p1err);
736 if ( currDraw%4 == 0 ){
737 currName = Form("vtx_%s_PtBin%i",fThetaAbsKeys->At(itheta)->GetName(), ibinpt);
738 can = new TCanvas(currName.Data(),currName.Data(),igroup1*xshift,igroup2*yshift,600,600);
741 can->cd( currDraw%4 + 1 );
744 fitFunc->DrawCopy("same");
747 SetSparseRange(gridSparse, kHvarMotherType, "", firstSrcBin, lastSrcBin, "USEBIN");
748 currName = Form("recoPt_%s",fThetaAbsKeys->At(itheta)->GetName());
749 can = new TCanvas(currName.Data(),currName.Data(),(igroup1+1)*xshift,igroup2*yshift,600,600);
750 TLegend* leg = new TLegend(0.6, 0.6, 0.8, 0.8);
752 for ( Int_t ireco=0; ireco<kNrecoHistos-1; ++ireco ) {
753 if ( recoHisto[ireco]->GetEntries() == 0. ) continue;
754 TH1* ratio = (TH1*)recoHisto[ireco]->Clone(Form("%s_ratio", recoHisto[ireco]->GetName()));
755 ratio->Divide(recoHisto[kRecoAll]);
756 ratio->SetLineColor(srcColors[ireco]);
757 ratio->SetMarkerColor(srcColors[ireco]);
758 ratio->SetMarkerStyle(20+ireco);
759 ratio->GetYaxis()->SetTitle("fraction of total");
760 ratio->Draw(drawOpt.Data());
761 leg->AddEntry(ratio,baseRecoName[ireco].Data(), "lp");
765 } // loop on theta abs