1 /**************************************************************************
2 * Copyright(c) 1998-2009, 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 **************************************************************************/
17 // Base class for DStar - Hadron Correlations Analysis
19 //-----------------------------------------------------------------------
22 // Author S.Bjelogrlic
24 // sandro.bjelogrlic@cern.ch
26 //-----------------------------------------------------------------------
30 //#include <TDatabasePDG.h>
31 #include <TParticle.h>
36 #include "AliAnalysisTaskDStarCorrelations.h"
37 #include "AliRDHFCutsDStartoKpipi.h"
38 #include "AliHFAssociatedTrackCuts.h"
39 #include "AliAODRecoDecay.h"
40 #include "AliAODRecoCascadeHF.h"
41 #include "AliAODRecoDecayHF2Prong.h"
42 #include "AliAODPidHF.h"
43 #include "AliVParticle.h"
44 #include "AliAnalysisManager.h"
45 #include "AliAODInputHandler.h"
46 #include "AliAODHandler.h"
47 #include "AliESDtrack.h"
48 #include "AliAODMCParticle.h"
49 #include "AliNormalizationCounter.h"
50 #include "AliReducedParticle.h"
51 #include "AliHFCorrelator.h"
57 ClassImp(AliAnalysisTaskDStarCorrelations)
60 //__________________________________________________________________________
61 AliAnalysisTaskDStarCorrelations::AliAnalysisTaskDStarCorrelations() :
80 // default constructor
83 //__________________________________________________________________________
84 AliAnalysisTaskDStarCorrelations::AliAnalysisTaskDStarCorrelations(const Char_t* name,AliRDHFCutsDStartoKpipi* cuts, AliHFAssociatedTrackCuts *AsscCuts) :
85 AliAnalysisTaskSE(name),
104 Info("AliAnalysisTaskDStarCorrelations","Calling Constructor");
105 DefineInput(0, TChain::Class());
106 DefineOutput(1,TList::Class()); // histos from data
107 DefineOutput(2,AliRDHFCutsDStartoKpipi::Class()); // my cuts
108 DefineOutput(3,AliNormalizationCounter::Class()); // normalization
109 DefineOutput(4,AliHFAssociatedTrackCuts::Class()); // my cuts
112 //__________________________________________________________________________
114 AliAnalysisTaskDStarCorrelations::~AliAnalysisTaskDStarCorrelations() {
119 Info("AliAnalysisTaskDStarCorrelations","Calling Destructor");
121 if(fhandler) {delete fhandler; fhandler = 0;}
122 //if(fPoolMgr) {delete fPoolMgr; fPoolMgr = 0;}
123 if(fmcArray) {delete fmcArray; fmcArray = 0;}
124 if(fCounter) {delete fCounter; fCounter = 0;}
125 if(fCorrelator) {delete fCorrelator; fCorrelator = 0;}
126 if(fOutput) {delete fOutput; fOutput = 0;}
127 if(fCuts) {delete fCuts; fCuts = 0;}
128 if(fCuts2) {delete fCuts2; fCuts2=0;}
132 //___________________________________________________________
133 void AliAnalysisTaskDStarCorrelations::Init(){
137 if(fDebug > 1) printf("AliAnalysisTaskDStarCorrelations::Init() \n");
139 AliRDHFCutsDStartoKpipi* copyfCuts=new AliRDHFCutsDStartoKpipi(*fCuts);
145 PostData(2,copyfCuts);
147 // Post the hadron cuts
156 //_________________________________________________
157 void AliAnalysisTaskDStarCorrelations::UserCreateOutputObjects(){
158 Info("UserCreateOutputObjects","CreateOutputObjects of task %s\n", GetName());
162 fOutput = new TList();
166 DefineHistoForAnalysis();
167 fCounter = new AliNormalizationCounter(Form("%s",GetOutputSlot(3)->GetContainer()->GetName()));
170 Double_t Pi = TMath::Pi();
171 fCorrelator = new AliHFCorrelator("Correlator",fCuts2,fSystem); // fCuts2 is the hadron cut object, fSystem to switch between pp or PbPb
172 fCorrelator->SetDeltaPhiInterval((-0.5-1./32)*Pi,(1.5-1./32)*Pi); // set correct phi interval
173 fCorrelator->SetEventMixing(fmixing); //set kFALSE/kTRUE for mixing Off/On
174 fCorrelator->SetAssociatedParticleType(fselect); // set 1/2/3 for hadron/kaons/kzeros
175 fCorrelator->SetApplyDisplacementCut(fDisplacement); //set kFALSE/kTRUE for using the displacement cut
176 fCorrelator->SetUseMC(fmontecarlo);
177 Bool_t pooldef = fCorrelator->DefineEventPool();
179 if(!pooldef) AliInfo("Warning:: Event pool not defined properly");
182 PostData(1,fOutput); // set the outputs
183 PostData(3,fCounter); // set the outputs
185 //_________________________________________________
186 void AliAnalysisTaskDStarCorrelations::UserExec(Option_t *){
190 std::cout << " " << std::endl;
191 std::cout << "=================================================================================" << std::endl;
193 if(fselect==1) std::cout << "TASK::Correlation with hadrons on SE "<< std::endl;
194 if(fselect==2) std::cout << "TASK::Correlation with kaons on SE "<< std::endl;
195 if(fselect==3) std::cout << "TASK::Correlation with kzeros on SE "<< std::endl;
198 if(fselect==1) std::cout << "TASK::Correlation with hadrons on ME "<< std::endl;
199 if(fselect==2) std::cout << "TASK::Correlation with kaons on ME "<< std::endl;
200 if(fselect==3) std::cout << "TASK::Correlation with kzeros on ME "<< std::endl;
203 Error("UserExec","NO EVENT FOUND!");
207 AliAODEvent* aodEvent = dynamic_cast<AliAODEvent*>(fInputEvent);
209 AliError("AOD event not found!");
213 fCorrelator->SetAODEvent(aodEvent); // set the event to be processed
215 fEvents++; // event counter
216 ((TH1D*)fOutput->FindObject("NofEvents"))->Fill(0);
217 fmcArray = dynamic_cast<TClonesArray*>(aodEvent->FindListObject(AliAODMCParticle::StdBranchName()));
219 if(fmontecarlo && !fmcArray){
220 AliError("Array of MC particles not found");
223 Bool_t isEvSel=fCuts->IsEventSelected(aodEvent);
225 ((TH1D*)fOutput->FindObject("NofEvents"))->Fill(1);
227 Bool_t correlatorON = fCorrelator->Initialize(); //define the pool for mixing
229 AliInfo("AliHFCorrelator didn't initialize the pool correctly or processed a bad event");
232 ((TH1D*)fOutput->FindObject("NofEvents"))->Fill(2);
234 if(fmontecarlo) fCorrelator->SetMCArray(fmcArray);
236 TClonesArray *arrayDStartoD0pi=0;
239 if(!aodEvent && AODEvent() && IsStandardAOD()) {
240 // In case there is an AOD handler writing a standard AOD, use the AOD
241 // event in memory rather than the input (ESD) event.
242 aodEvent = dynamic_cast<AliAODEvent*> (AODEvent());
243 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
244 // have to taken from the AOD event hold by the AliAODExtension
245 AliAODHandler* aodHandler = (AliAODHandler*)
246 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
247 if(aodHandler->GetExtensions()) {
248 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
249 AliAODEvent *aodFromExt = ext->GetAOD();
250 arrayDStartoD0pi=(TClonesArray*)aodFromExt->GetList()->FindObject("Dstar");
253 arrayDStartoD0pi=(TClonesArray*)aodEvent->GetList()->FindObject("Dstar");
256 if(!aodEvent->GetPrimaryVertex() || TMath::Abs(aodEvent->GetMagneticField())<0.001) return;
258 // initialize variables you will need for the D*
263 Bool_t isInPeak, isInSideBand, isDStarMCtag;
264 Double_t invMassDZero;
265 Double_t deltainvMDStar;
268 Double_t mPDGD0=1.8648;//TDatabasePDG::Instance()->GetParticle(421)->Mass();
269 Double_t mPDGDstar=2.01022;//TDatabasePDG::Instance()->GetParticle(413)->Mass();
272 //MC tagging for DStar
273 //D* and D0 prongs needed to MatchToMC method
274 Int_t pdgDgDStartoD0pi[2]={421,211};
275 Int_t pdgDgD0toKpi[2]={321,211};
278 //loop on D* candidates
279 for (Int_t iDStartoD0pi = 0; iDStartoD0pi<arrayDStartoD0pi->GetEntriesFast(); iDStartoD0pi++) {
281 isInSideBand = kFALSE;
282 isDStarMCtag = kFALSE;
286 invMassDZero = - 999;
287 deltainvMDStar = -998;
289 AliAODRecoCascadeHF* dstarD0pi = (AliAODRecoCascadeHF*)arrayDStartoD0pi->At(iDStartoD0pi);
290 if(!dstarD0pi->GetSecondaryVtx()) continue;
291 AliAODRecoDecayHF2Prong* theD0particle = (AliAODRecoDecayHF2Prong*)dstarD0pi->Get2Prong();
292 if (!theD0particle) continue;
295 // track quality cuts
296 Int_t isTkSelected = fCuts->IsSelected(dstarD0pi,AliRDHFCuts::kTracks); // quality cuts on tracks
297 // region of interest + topological cuts + PID
298 Int_t isSelected=fCuts->IsSelected(dstarD0pi,AliRDHFCuts::kCandidate); //selected
300 if(!isTkSelected) continue;
301 if(!isSelected) continue;
302 if(!fCuts->IsInFiducialAcceptance(dstarD0pi->Pt(),dstarD0pi->YDstar())) continue;
303 Int_t mcLabelDStar = -999;
305 // find associated MC particle for D* ->D0toKpi
306 mcLabelDStar = dstarD0pi->MatchToMC(413,421,pdgDgDStartoD0pi,pdgDgD0toKpi,fmcArray,kFALSE);
307 if(mcLabelDStar>=0) isDStarMCtag = kTRUE;
310 ptDStar = dstarD0pi->Pt();
311 phiDStar = dstarD0pi->Phi();
312 etaDStar = dstarD0pi->Eta();
314 phiDStar = fCorrelator->SetCorrectPhiRange(phiDStar); // set the phi of the D meson in the correct range
316 Int_t ptbin=fCuts->PtBin(dstarD0pi->Pt());
318 Double_t dmDStarWindow =0.0019;// 0.0019 = 3 sigma
319 Double_t mD0Window=0.074;
322 if (ptbin==1) mD0Window = 0.026; //0.5-1
323 if (ptbin==2) mD0Window = 0.022; //1-2
324 if (ptbin==3) mD0Window = 0.024; //2-3
325 if (ptbin==4) mD0Window = 0.032;
326 if (ptbin==5) mD0Window = 0.032;
327 if (ptbin==6) mD0Window = 0.036;
328 if (ptbin==7) mD0Window = 0.036;
329 if (ptbin==8) mD0Window = 0.036;
330 if (ptbin==9) mD0Window = 0.058;
331 if (ptbin==10) mD0Window = 0.058;
332 if (ptbin>10) mD0Window = 0.074;
335 if (ptbin==0) mD0Window = 0.032; //1-1
336 if (ptbin==1) mD0Window = 0.032; //2-3
337 if (ptbin==2) mD0Window = 0.032; //3-4
338 if (ptbin==3) mD0Window = 0.032; //4-5
339 if (ptbin==4) mD0Window = 0.036; //5-6
340 if (ptbin==5) mD0Window = 0.036; //6-8
341 if (ptbin==6) mD0Window = 0.055; //8-12
342 if (ptbin==7) mD0Window = 0.074; //12-16
343 if (ptbin==8) mD0Window = 0.074; //16-24
344 if (ptbin==9) mD0Window = 0.074; //24-35
347 invMassDZero = dstarD0pi->InvMassD0();
348 ((TH2F*)fOutput->FindObject("D0InvMass"))->Fill(ptDStar,invMassDZero);
350 deltainvMDStar = dstarD0pi->DeltaInvMass();
354 if (TMath::Abs(invMassDZero-mPDGD0)<mD0Window){
356 ((TH2F*)fOutput->FindObject("DeltaInvMass"))->Fill(ptDStar,deltainvMDStar);
357 if(TMath::Abs(deltainvMDStar-(mPDGDstar-mPDGD0))<dmDStarWindow){ // is in DStar peak region?
359 ((TH1F*)fOutput->FindObject("RecoPtDStar"))->Fill(ptDStar);
361 ((TH2F*)fOutput->FindObject("PhiEtaTrigger"))->Fill(phiDStar,etaDStar);
363 }// end if good candidates
366 if (TMath::Abs(invMassDZero-mPDGD0)>1.3*mD0Window && TMath::Abs(invMassDZero-mPDGD0)<4.*mD0Window ){
367 ((TH2F*)fOutput->FindObject("bkgDeltaInvMass"))->Fill(ptDStar,deltainvMDStar);
368 ((TH2F*)fOutput->FindObject("D0InvMassinSB"))->Fill(ptDStar,invMassDZero);
370 if(TMath::Abs(deltainvMDStar-(mPDGDstar-mPDGD0))<dmDStarWindow){ // is in DStar peak region?
371 ((TH1F*)fOutput->FindObject("RecoPtBkg"))->Fill(ptDStar);
372 isInSideBand = kTRUE;
373 ((TH2F*)fOutput->FindObject("PhiEtaSideBand"))->Fill(phiDStar,etaDStar);
377 // getting the number of triggers in the MCtag D* case
380 if(fmontecarlo && isDStarMCtag) ((TH1F*)fOutput->FindObject("MCtagPtDStar"))->Fill(ptDStar);
383 if(!isInPeak && !isInSideBand) continue; // skip if it is not side band or peak event - SAVE CPU TIME
386 fCorrelator->SetTriggerParticleProperties(ptDStar,phiDStar,etaDStar); // pass to the object the necessary trigger part parameters
388 Short_t daughtercharge = ((AliAODTrack*)theD0particle->GetDaughter(0))->Charge();
389 fCorrelator->SetTriggerParticleDaughterCharge(daughtercharge);
392 Int_t trackiddaugh0 = ((AliAODTrack*)theD0particle->GetDaughter(0))->GetID();
393 Int_t trackiddaugh1 = ((AliAODTrack*)theD0particle->GetDaughter(1))->GetID();
394 Int_t trackidsoftPi = ((AliAODTrack*)dstarD0pi->GetBachelor())->GetID();
396 Bool_t execPool = fCorrelator->ProcessEventPool();
397 if(fmixing && !execPool) {
398 AliInfo("Mixed event analysis: pool is not ready");
402 Int_t NofEventsinPool = 1;
403 if(fmixing) NofEventsinPool = fCorrelator->GetNofEventsInPool();
405 for (Int_t jMix =0; jMix < NofEventsinPool; jMix++){// loop on events in the pool; if it is SE analysis, stops at one
407 Bool_t analyzetracks = fCorrelator->ProcessAssociatedTracks(jMix);
410 AliInfo("AliHFCorrelator::Cannot process the track array");
415 Int_t NofTracks = fCorrelator->GetNofTracks();
417 for(Int_t iTrack = 0; iTrack<NofTracks; iTrack++){ // looping on track candidates
419 Bool_t runcorrelation = fCorrelator->Correlate(iTrack);
420 if(!runcorrelation) continue;
422 Double_t DeltaPhi = fCorrelator->GetDeltaPhi();
423 Double_t DeltaEta = fCorrelator->GetDeltaEta();
425 AliReducedParticle * hadron = fCorrelator->GetAssociatedParticle();
427 Double_t ptHad = hadron->Pt();
428 Double_t phiHad = hadron->Phi();
429 Double_t etaHad = hadron->Eta();
430 Double_t label = hadron->GetLabel();
431 Double_t trackid = hadron->GetID();
433 phiHad = fCorrelator->SetCorrectPhiRange(phiHad);
435 if(!fmixing){ // skip D* Daughetrs
436 if(trackid == trackiddaugh0) continue;
437 if(trackid == trackiddaugh1) continue;
438 if(trackid == trackidsoftPi) continue;
441 // from here on it is up to the user to decide what object to fill
443 if(fmontecarlo && isDStarMCtag){ // check correlations of MC tagged DStars in MonteCarlo
445 Int_t PartSource = fCuts2->IsMCpartFromHF(label,fmcArray); // check source of associated particle (hadron/kaon/K0)
447 FillMCTagCorrelations(ptDStar,DeltaPhi,DeltaEta,ptHad,PartSource);
453 if(fselect==1) ((TH3D*)fOutput->FindObject("DPhiDStarHadron"))->Fill(DeltaPhi,ptDStar,DeltaEta);
454 if(fselect==2) ((TH3D*)fOutput->FindObject("DPhiDStarKaon"))->Fill(DeltaPhi,ptDStar,DeltaEta);
455 if(fselect==3) ((TH3D*)fOutput->FindObject("DPhiDStarKZero"))->Fill(DeltaPhi,ptDStar,DeltaEta);
456 ((TH2F*)fOutput->FindObject("PhiEtaPart"))->Fill(phiHad,etaHad);
457 //counterPeak++; // count tracks per peak per event
463 if(fselect==1) ((TH3D*)fOutput->FindObject("bkgDPhiDStarHadron"))->Fill(DeltaPhi,ptDStar,DeltaEta);
464 if(fselect==2) ((TH3D*)fOutput->FindObject("bkgDPhiDStarKaon"))->Fill(DeltaPhi,ptDStar,DeltaEta);
465 if(fselect==3) ((TH3D*)fOutput->FindObject("bkgDPhiDStarKZero"))->Fill(DeltaPhi,ptDStar,DeltaEta);
473 } // end loop on track candidates
474 } // end loop on events in the pool
476 }// end loop on D* candidates
477 Bool_t updated = fCorrelator->PoolUpdate();
478 if(!updated) AliInfo("Pool was not updated");
480 //std::cout << ">>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> END OF THE EVENT <<<<<<<<<<<<<<<<<<<<<<<<<<<<<" << std::endl;
486 //________________________________________ terminate ___________________________
487 void AliAnalysisTaskDStarCorrelations::Terminate(Option_t*)
489 // The Terminate() function is the last function to be called during
490 // a query. It always runs on the client, it can be used to present
491 // the results graphically or save the results to file.
493 AliAnalysisTaskSE::Terminate();
495 fOutput = dynamic_cast<TList*> (GetOutputData(1));
497 printf("ERROR: fOutput not available\n");
505 //_____________________________________________________
506 void AliAnalysisTaskDStarCorrelations::DefineHistoForAnalysis(){
508 Double_t Pi = TMath::Pi();
509 Int_t nbinscorr = 32;
510 Double_t lowcorrbin = -0.5*Pi - Pi/32; // shift the bin by half the width so that at 0 is it the bin center
511 Double_t upcorrbin = 1.5*Pi - Pi/32;
513 // ========================= histograms for both Data and MonteCarlo
516 TH1D * NofEvents = new TH1D("NofEvents","NofEvents",3,0,3);
517 fOutput->Add(NofEvents);
519 TH2F *D0InvMass = new TH2F("D0InvMass","K#pi invariant mass distribution",300,0,30,1500,0.5,3.5);
520 fOutput->Add(D0InvMass);
522 TH2F *D0InvMassinSB = new TH2F("D0InvMassinSB","K#pi invariant mass distribution in sb",300,0,30,1500,0.5,3.5);
523 fOutput->Add(D0InvMassinSB);
525 TH2F *DeltaInvMass = new TH2F("DeltaInvMass","K#pi#pi - K#pi invariant mass distribution",300,0,30,750,0.1,0.2);
526 fOutput->Add(DeltaInvMass);
528 TH2F *bkgDeltaInvMass = new TH2F("bkgDeltaInvMass","K#pi#pi - K#pi invariant mass distribution",300,0,30,750,0.1,0.2);
529 fOutput->Add(bkgDeltaInvMass);
531 TH1F *RecoPtDStar = new TH1F("RecoPtDStar","RECO DStar pt distribution",30,0,30);
532 fOutput->Add(RecoPtDStar);
534 TH1F *RecoPtBkg = new TH1F("RecoPtBkg","RECO pt distribution side bands",30,0,30);
535 fOutput->Add(RecoPtBkg);
537 TH1F *MCtagPtDStar = new TH1F("MCtagPtDStar","RECO pt of MCtagged DStars side bands",30,0,30);
538 fOutput->Add(MCtagPtDStar);
540 TH2F *KZeroSpectra = new TH2F("KZeroSpectra","Spectra of K0s",500,0.3,0.8,250,0,25);
541 if(fselect==3) fOutput->Add(KZeroSpectra);
543 TH2F *KZeroSpectraifHF = new TH2F("KZeroSpectraifHF","Spectra of K0s in association with a D*",500,0.3,0.8,250,0,25);
544 if(fselect==3) fOutput->Add(KZeroSpectraifHF);
546 TH1D * NofTracksInPeak = new TH1D("NofTracksInPeak","NofTracksInPeak",500,0.5,500.5);
547 fOutput->Add(NofTracksInPeak);
549 TH1D * NofTracksInSB = new TH1D("NofTracksInSB","NofTracksInSB",500,0.5,500.5);
550 fOutput->Add(NofTracksInSB);
552 TH2I * EventMixingCheck = new TH2I("EventMixingCheck","EventMixingCheck",5,-0.5,4.5,7,-0.5,6.5);
553 if(fmixing) fOutput->Add(EventMixingCheck);
557 TH1F * MCSources = new TH1F("MCSources","Origin of associated particles in MC", 10, -0.5, 9.5);
558 MCSources->GetXaxis()->SetBinLabel(1,"All ");
559 MCSources->GetXaxis()->SetBinLabel(2," from Heavy flavour");
560 MCSources->GetXaxis()->SetBinLabel(3," from c->D");
561 MCSources->GetXaxis()->SetBinLabel(4," from b->D");
562 MCSources->GetXaxis()->SetBinLabel(5," from b->B");
563 if(fmontecarlo) fOutput->Add(MCSources);
565 TH2F * PhiEtaTrigger = new TH2F("PhiEtaTrigger","#phi distribution of the trigger particle",36,-0.5*Pi,1.5*Pi,18,-0.9,0.9);
566 fOutput->Add(PhiEtaTrigger);
568 TH2F * PhiEtaSideBand = new TH2F("PhiEtaSideBand","#phi distribution of the sideband particle",36,-0.5*Pi,1.5*Pi,18,-0.9,0.9);
569 fOutput->Add(PhiEtaSideBand);
571 TH2F * PhiEtaPart = new TH2F("PhiEtaPart","#phi distribution of the associated particle",36,-0.5*Pi,1.5*Pi,18,-0.9,0.9);
572 fOutput->Add(PhiEtaPart);
575 //correlations histograms
576 TString histoname1 = "DPhiDStar";
577 if(fselect==1) histoname1 += "Hadron";
578 if(fselect==2) histoname1 += "Kaon";
579 if(fselect==3) histoname1 += "KZero";
582 TH3D * DPhiDStar = new TH3D(histoname1.Data(),histoname1.Data(),nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
584 TH3D * DPhiDStarKZero1 = new TH3D("DPhiDStarKZero1","DPhiDStarKZero1",nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
586 //side band background histograms
587 TString histoname2 = "bkg";
588 histoname2 += histoname1;
589 TH3D * bkgDPhiDStar = new TH3D(histoname2.Data(),histoname2.Data(),nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
590 TH3D * bkgDPhiDStarKZero1 = new TH3D("bkgDPhiDStarKZero1","bkgDPhiDStarKZero1",nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
593 fOutput->Add(DPhiDStar);
595 if(fselect==3){fOutput->Add(DPhiDStarKZero1);}
597 fOutput->Add(bkgDPhiDStar);
599 if(fselect==3){fOutput->Add(bkgDPhiDStarKZero1);}
602 // ========================= histos for analysis on MC
603 TString histoname3 = "MCTag";
604 histoname3 += histoname1;
605 TH3D * MCTagDPhiDStar = new TH3D(histoname3.Data(),histoname3.Data(),nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
607 TString histoname44 = "CharmDOrigin";
608 histoname44 += histoname1;
611 TH3D * CharmDOriginDPhiDStar = new TH3D(histoname44.Data(),histoname44.Data(),nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
614 TString histoname54 = "BeautyDOrigin";
615 histoname54 += histoname1;
617 TH3D * BeautyDOriginDPhiDStar = new TH3D(histoname54.Data(),histoname54.Data(),nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
619 TString histoname55 = "BeautyBOrigin";
620 histoname55 += histoname1;
622 TH3D * BeautyBOriginDPhiDStar = new TH3D(histoname55.Data(),histoname55.Data(),nbinscorr,lowcorrbin,upcorrbin,30,0,30,19,-0.95,0.95);
626 fOutput->Add(MCTagDPhiDStar);
627 fOutput->Add(CharmDOriginDPhiDStar);
628 fOutput->Add(BeautyDOriginDPhiDStar);
629 fOutput->Add(BeautyBOriginDPhiDStar);
639 //____________________________ Function for MC correlations ___________________________________________________
640 void AliAnalysisTaskDStarCorrelations::FillMCTagCorrelations(Double_t ptTrig, Double_t DelPhi, Double_t DelEta, Double_t ptTrack, Int_t mcSource){
643 if(fselect==1) ((TH3D*)fOutput->FindObject("MCTagDPhiDStarHadron"))->Fill(DelPhi,ptTrig,DelEta);
644 if(fselect==2 && ptTrack <1.5) ((TH3D*)fOutput->FindObject("MCTagDPhiDStarKaon"))->Fill(DelPhi,ptTrig,DelEta);
645 if(fselect==3) ((TH3D*)fOutput->FindObject("MCTagDPhiDStarKZero"))->Fill(DelPhi,ptTrig,DelEta);
649 ((TH1F*)fOutput->FindObject("MCSources"))->Fill(0);
651 if (mcSource==44){ // is from charm ->D
652 if(fselect==1) ((TH3D*)fOutput->FindObject("CharmDOriginDPhiDStarHadronMC"))->Fill(DelPhi,ptTrig,DelEta);
653 if(fselect==2 && ptTrack <1.5) ((TH3D*)fOutput->FindObject("CharmDOriginDPhiDStarKaonMC"))->Fill(DelPhi,ptTrig,DelEta);
654 if(fselect==3) ((TH3D*)fOutput->FindObject("CharmDOriginDPhiDStarKZeroMC"))->Fill(DelPhi,ptTrig,DelEta);
657 ((TH1F*)fOutput->FindObject("MCSources"))->Fill(1);
658 ((TH1F*)fOutput->FindObject("MCSources"))->Fill(2);
661 if (mcSource==54){ // is from beauty -> D
662 if(fselect==1) ((TH3D*)fOutput->FindObject("BeautyDOriginDPhiDStarHadronMC"))->Fill(DelPhi,ptTrig,DelEta);
663 if(fselect==2 && ptTrack <1.5) ((TH3D*)fOutput->FindObject("BeautyDOriginDPhiDStarKaonMC"))->Fill(DelPhi,ptTrig,DelEta);
664 if(fselect==3) ((TH3D*)fOutput->FindObject("BeautyDOriginDPhiDStarKZeroMC"))->Fill(DelPhi,ptTrig,DelEta);
665 if(fselect==3) ((TH1F*)fOutput->FindObject("MCSources"))->Fill(1);
666 if(fselect==3) ((TH1F*)fOutput->FindObject("MCSources"))->Fill(3);
669 if (mcSource==55){ // is from beauty ->B
670 if(fselect==1) ((TH3D*)fOutput->FindObject("BeautyBOriginDPhiDStarHadronMC"))->Fill(DelPhi,ptTrig,DelEta);
671 if(fselect==2 && ptTrack <1.5) ((TH3D*)fOutput->FindObject("BeautyBOriginDPhiDStarKaonMC"))->Fill(DelPhi,ptTrig,DelEta);
672 if(fselect==3) ((TH3D*)fOutput->FindObject("BeautyOriginBDPhiDStarKZeroMC"))->Fill(DelPhi,ptTrig,DelEta);
673 if(fselect==3) ((TH1F*)fOutput->FindObject("MCSources"))->Fill(1);
674 if(fselect==3) ((TH1F*)fOutput->FindObject("MCSources"))->Fill(4);