ROOT logo
/**************************************************************************
* Copyright(c) 1998-2008,ALICE Experiment at CERN,All rights reserved.    *
*                                                                         *
* Author: The ALICE Off-line Project.                                     *
* Contributors are mentioned in the code where appropriate.               *
*                                                                         *
* Permission to use,copy,modify and distribute this software and its      *
* documentation strictly for non-commercial purposes is hereby granted    *
* without fee,provided that the above copyright notice appears in all     *
* copies and that both the copyright notice and this permission notice    *
* appear in the supporting documentation. The authors make no claims      *
* about the suitability of this software for any purpose. It is           *
* provided "as is" without express or implied warranty.                   *
**************************************************************************/

//==============================================================================
// FlowD2H main task:
// >> Select candidates and passes flowevents to the daughter tasks.
// >> The POIcuts are polymorphic based on the AliRDHFCuts class allowing the 
//    use of all charmed candidates reconstructed in the central barrel.
// Author: Carlos Perez (cperez@cern.ch)
//==============================================================================

/* $Id$ */

#include "TChain.h"

#include "AliAnalysisTaskSE.h"
#include "AliAnalysisManager.h"

#include "AliCentrality.h"

#include "AliAODEvent.h"
#include "AliAODInputHandler.h"
#include "AliAODTrack.h"

#include "AliRDHFCuts.h"

#include "AliRDHFCutsD0toKpi.h"
#include "AliRDHFCutsDStartoKpipi.h"
#include "AliRDHFCutsDplustoKpipi.h"

#include "AliAODRecoDecayHF2Prong.h"
#include "AliAODRecoDecayHF3Prong.h"
#include "AliAODRecoCascadeHF.h"

#include "TMath.h"

#include "AliFlowCommonConstants.h"
#include "AliFlowEvent.h"
#include "AliFlowTrackCuts.h"
#include "AliFlowCandidateTrack.h"

#include "TObjArray.h"
#include "TList.h"
#include "TH1D.h"

#include "AliAnalysisTaskFlowD2H.h"

ClassImp(AliAnalysisTaskFlowD2H)

//_____________________________________________________________________________
AliAnalysisTaskFlowD2H::AliAnalysisTaskFlowD2H() :
  AliAnalysisTaskSE(), fTPCEvent(NULL), fVZEEvent(NULL), 
  fCutsTPC(NULL), fCutsVZE(NULL), fNoPOIs(NULL), fCutsPOI(NULL),
  fSource(0), fDebugV2(kFALSE), fSwap(kFALSE), fMassBins(0), fMinMass(0.),
  fMaxMass(0.), fPtBinWidth(0), fHList(NULL), fEvent(NULL), 
  fCC(NULL), fRFPMTPC(NULL), fRFPPhiTPC(NULL), fCandidates(NULL)
{
// Default constructor
}

//_____________________________________________________________________________
AliAnalysisTaskFlowD2H::AliAnalysisTaskFlowD2H(const char *name,	
					       AliFlowTrackCuts *cutsTPC,
					       AliFlowTrackCuts *cutsVZE,
					       AliRDHFCuts *cutsPOIs,
					       Int_t specie) :
  AliAnalysisTaskSE(name), fTPCEvent(NULL), fVZEEvent(NULL), 
  fCutsTPC(cutsTPC), fCutsVZE(cutsVZE), fNoPOIs(NULL), fCutsPOI(cutsPOIs),
  fSource(specie), fDebugV2(kFALSE), fSwap(kFALSE), fMassBins(0), fMinMass(0.),
  fMaxMass(0.), fPtBinWidth(0), fHList(NULL), fEvent(NULL),
  fCC(NULL), fRFPMTPC(NULL), fRFPPhiTPC(NULL), fCandidates(NULL)
{
// Standard constructor
  DefineInput( 0,TChain::Class());
  DefineOutput(1,TList::Class());
  DefineOutput(2,AliFlowEventSimple::Class());
  DefineOutput(3,AliFlowEventSimple::Class());
}

//_____________________________________________________________________________
AliAnalysisTaskFlowD2H::~AliAnalysisTaskFlowD2H(){
  // delete objects
  if(fTPCEvent) delete fTPCEvent;
  if(fVZEEvent) delete fVZEEvent;
  if(fCutsTPC) delete fCutsTPC;
  if(fCutsVZE) delete fCutsVZE;
  if(fNoPOIs) delete fNoPOIs;
  if(fCutsPOI) delete fCutsPOI;
  if(fHList) delete fHList;
  if(fCandidates) delete fCandidates;
}

//_____________________________________________________________________________
void AliAnalysisTaskFlowD2H::UserCreateOutputObjects(){
  // Define output objects + parameters
  if(fDebugV2)
    printf("DEBUGGER: Creating output\n");
  fHList = new TList();
  fHList->SetOwner();
  fEvent = new TH1D("Events","Events",3,0,3);
  fEvent->GetXaxis()->SetBinLabel(1,"REACHED");
  fEvent->GetXaxis()->SetBinLabel(2,"SELECTED");
  fEvent->GetXaxis()->SetBinLabel(3,"DELTA AOD REACHED");
  fHList->Add( fEvent );
  fCC = new TH1D("CentralityClass","Centrality Class",50,0,100);
  fHList->Add( fCC );
  fRFPMTPC = new TH1D("RFPMultiplicityTPC","RFP Multiplicity TPC",300,0,3000);
  fHList->Add( fRFPMTPC );
  fRFPPhiTPC = new TH1D("RFPPhiTPC","RFP Phi TPC",180,0,TMath::TwoPi());
  fHList->Add( fRFPPhiTPC );

  fCandidates = new TObjArray(100);
  fCandidates->SetOwner();

  AliFlowCommonConstants* cc = AliFlowCommonConstants::GetMaster();
  cc->SetNbinsMult(500);
  cc->SetNbinsPt(24/fPtBinWidth);
  cc->SetNbinsPhi(360);
  cc->SetNbinsEta(90);
  cc->SetNbinsQ(100);
  cc->SetNbinsMass( fMassBins );
  cc->SetMultMin(0);
  cc->SetMultMax(500);
  cc->SetPtMin(0);
  cc->SetPtMax(24);
  cc->SetPhiMin(0);
  cc->SetPhiMax(TMath::TwoPi());
  cc->SetEtaMin(-3.9);
  cc->SetEtaMax(+5.1);
  cc->SetQMin(0);
  cc->SetQMax(1);
  cc->SetMassMin( fMinMass );
  cc->SetMassMax( fMaxMass );

  fTPCEvent = new AliFlowEvent(3000);
  fVZEEvent = new AliFlowEvent(170);

  fNoPOIs = new AliFlowTrackCuts( "noPOIs" );
  fNoPOIs->SetParamType(AliFlowTrackCuts::kGlobal);
  fNoPOIs->SetPtRange(+1,-1);

  PostData(1,fHList);
  PostData(2,fTPCEvent);
  PostData(3,fVZEEvent);
}
//_____________________________________________________________________________
void AliAnalysisTaskFlowD2H::UserExec(Option_t *)
{
  // Do analysis + fIll histograms
  AliAODEvent *fAOD = dynamic_cast<AliAODEvent*>(InputEvent());

  if(!fAOD) return;
  fEvent->Fill( 0 );

  // floweventcuts::isselected() and alirdhfcuts::iseventselected() cut on the same 
  // values in the same way BUT the latter also loads the PIDresponse object from the 
  // event header!!! 
  //  if(!fEventCuts->IsSelected(fAOD)) return;
  if(!fCutsPOI->IsEventSelected(fAOD)) return;
  fEvent->Fill( 1 );

  fCC->Fill( fCutsPOI->GetCentrality(fAOD) );

  fCutsTPC->SetEvent( fAOD, MCEvent() );
  fCutsVZE->SetEvent( fAOD, MCEvent() );
  fNoPOIs->SetEvent( fAOD, MCEvent() );
  fTPCEvent->Fill( fCutsTPC, fNoPOIs );
  fVZEEvent->Fill( fCutsVZE, fNoPOIs );

  Int_t rfps = fTPCEvent->GetNumberOfRPs();
  fRFPMTPC->Fill( rfps );
  for(int iRPs=0; iRPs!=rfps; ++iRPs ) {
    AliFlowTrack *iRP = dynamic_cast<AliFlowTrack*>(fTPCEvent->GetTrack( iRPs ));
    if (!iRP) continue;
    fRFPPhiTPC->Fill( iRP->Phi() );
  }

  if (fDebugV2) printf("Event selected\n");
  fCandidates->SetLast(-1); // resets the array

  switch (fSource) {
    case (AliRDHFCuts::kD0toKpiCuts):
      FillD0toKpi(fAOD); break;
    case (AliRDHFCuts::kDstarCuts):
      FillDStartoKpipi(fAOD); break;
    case (AliRDHFCuts::kDplusCuts):
      FillDplustoKpipi(fAOD); break;
  }

  if(fDebugV2) printf("TPCevent %d | VZEevent %d\n", fTPCEvent->NumberOfTracks(), fVZEEvent->NumberOfTracks() );
  //inject candidates
  if (fDebugV2)  printf("I received %d candidates\n",fCandidates->GetEntriesFast());
  for(int iCand=0; iCand!=fCandidates->GetEntriesFast(); ++iCand ) {
    AliFlowCandidateTrack *cand = dynamic_cast<AliFlowCandidateTrack*>(fCandidates->At(iCand));
    if (!cand) continue;
    if (fDebugV2) printf(" >Checking at candidate %d with %d daughters: mass %f\n",iCand,cand->GetNDaughters(),cand->Mass());
    for(int iDau=0; iDau!=cand->GetNDaughters(); ++iDau) {
      if(fDebugV2) printf("  >Daughter %d with fID %d", iDau, cand->GetIDDaughter(iDau));
      for(int iRPs=0; iRPs!=fTPCEvent->NumberOfTracks(); ++iRPs ) {
	AliFlowTrack *iRP = dynamic_cast<AliFlowTrack*>(fTPCEvent->GetTrack( iRPs ));
	if (!iRP) continue;
	if( !iRP->InRPSelection() ) continue;
	if( cand->GetIDDaughter(iDau) == iRP->GetID() ) {
	  if(fDebugV2) printf(" was in RP set");
	  iRP->SetForRPSelection(kFALSE);
	  fTPCEvent->SetNumberOfRPs( fTPCEvent->GetNumberOfRPs() -1 );
	}
      }
      if(fDebugV2) printf("\n");
    }
    cand->SetForPOISelection(kTRUE);
    fTPCEvent->InsertTrack( ((AliFlowTrack*) cand) );
    fVZEEvent->InsertTrack( ((AliFlowTrack*) cand) );
    fTPCEvent->SetNumberOfPOIs(fTPCEvent->GetNumberOfPOIs()+1);
    fVZEEvent->SetNumberOfPOIs(fVZEEvent->GetNumberOfPOIs()+1);
  }
  if(fDebugV2) printf("TPCevent %d | VZEevent %d\n", fTPCEvent->NumberOfTracks(), fVZEEvent->NumberOfTracks() );

  PostData(1,fHList);
  PostData(2,fTPCEvent);
  PostData(3,fVZEEvent);
}
//______________________________________________________________________________
void AliAnalysisTaskFlowD2H::FillD0toKpi(const AliAODEvent *theAOD)
{
  // Fill D0->Kpi histos
  TList *listHF = (TList*) theAOD->GetList();
  if(!listHF) return;
  TClonesArray *listDzero = (TClonesArray*) listHF->FindObject("D0toKpi");
  if(!listDzero) return;
  int nEntries = listDzero->GetEntriesFast();
  if( !nEntries ) return;
  fEvent->Fill( 2 ); // EVERYTHING OKAY

  AliRDHFCutsD0toKpi *fCutsD0toKpi = (AliRDHFCutsD0toKpi*) fCutsPOI;
  if (fDebugV2) printf("  ᶫ%d candidates found. Looping...\n", nEntries);
  Int_t nIDs[2];
  for( int iEntry=0; iEntry!=nEntries; ++iEntry ) {
    AliAODRecoDecayHF2Prong *d0cand = (AliAODRecoDecayHF2Prong*) listDzero->UncheckedAt( iEntry );
    if( !d0cand ) continue;
    // APPLYING CUTS
    if( !d0cand->HasSelectionBit(AliRDHFCuts::kD0toKpiCuts) ) continue;
    if( !fCutsD0toKpi->IsInFiducialAcceptance(d0cand->Pt(),d0cand->Y(421)) )continue;
    Int_t topCut = fCutsD0toKpi->IsSelected( d0cand, AliRDHFCuts::kAll, NULL );
    if(!topCut) continue;
    Double_t massD0=topCut>1?d0cand->InvMassD0bar():d0cand->InvMassD0();
    if((fSwap)&&(topCut==3)) {
      massD0=d0cand->InvMassD0();
    }

    // TO HANDLE AUTOCORRELATIONS
    nIDs[0] = ( (AliAODTrack*)d0cand->GetDaughter(0) )->GetID();
    nIDs[1] = ( (AliAODTrack*)d0cand->GetDaughter(1) )->GetID();
    // ADDING TRACK
    MakeTrack(massD0, d0cand->Pt(), d0cand->Phi(), d0cand->Eta(), 2, nIDs);
    if(fDebugV2) printf("   ᶫInjecting D0 candidate\n");
  }
}
//______________________________________________________________________________
void AliAnalysisTaskFlowD2H::FillDStartoKpipi(const AliAODEvent *theAOD )
{
  // Fills D* to Kpipi
  TList *listHF = (TList*) theAOD->GetList();
  if(!listHF) return;
  TClonesArray *listDstar = (TClonesArray*) listHF->FindObject("Dstar");
  if(!listDstar) return;
  int nEntries = listDstar->GetEntriesFast();
  if( !nEntries ) return;
  fEvent->Fill( 2 ); // EVERYTHING OKAY

  AliRDHFCutsDStartoKpipi *fCutsDStartoKpipi = (AliRDHFCutsDStartoKpipi*) fCutsPOI;
  if (fDebugV2) printf("  ᶫ%d candidates found. Looping...\n", nEntries);
  Int_t nIDs[3];
  for( int iEntry=0; iEntry!=nEntries; ++iEntry ) {
    AliAODRecoCascadeHF *dst = (AliAODRecoCascadeHF*) listDstar->UncheckedAt( iEntry );
    if( !dst ) continue;
    AliAODRecoDecayHF2Prong *d0cand = (AliAODRecoDecayHF2Prong*)dst->Get2Prong();
    if(!d0cand) continue;
    if( !d0cand->GetDaughter(0) ) continue;
    if( !d0cand->GetDaughter(1) ) continue;
    if( !dst->GetBachelor() ) continue;
    // APPLYING CUTS
    if( !dst->HasSelectionBit(AliRDHFCuts::kDstarCuts) ) continue;
    if( !fCutsDStartoKpipi->IsInFiducialAcceptance(dst->Pt(),dst->YDstar()) )continue;
    Int_t topCut=0;
    topCut = fCutsDStartoKpipi->IsSelected( dst, AliRDHFCuts::kAll );
    if(!topCut) continue;
    Double_t massDS=dst->DeltaInvMass();
    // TO HANDLE AUTOCORRELATIONS
    nIDs[0] = ((AliAODTrack*)d0cand->GetDaughter(0))->GetID();
    nIDs[1] = ((AliAODTrack*)d0cand->GetDaughter(1))->GetID();
    nIDs[2] = ((AliAODTrack*)dst->GetBachelor() )->GetID();
    // ADDING TRACK
    MakeTrack(massDS, dst->Pt(), dst->Phi(), dst->Eta(), 3, nIDs);
    if(fDebugV2) printf("   ᶫInjecting DStar candidate %d\n",iEntry);
  }
}
//______________________________________________________________________________
void AliAnalysisTaskFlowD2H::FillDplustoKpipi(const AliAODEvent *theAOD )
{
  // Fill D+ to Kpipi
  TList *listHF = (TList*) theAOD->GetList();
  if(!listHF) return;
  TClonesArray *listDplus = (TClonesArray*) listHF->FindObject("Charm3Prong");
  if(!listDplus) return;
  int nEntries = listDplus->GetEntriesFast();
  if( !nEntries ) return;
  fEvent->Fill( 2 ); // EVERYTHING OKAY

  AliRDHFCutsDplustoKpipi *fCutsDStartoKpipi = (AliRDHFCutsDplustoKpipi*) fCutsPOI;
  if (fDebugV2) printf("  ᶫ%d candidates found. Looping...\n", nEntries);
  Int_t nIDs[3];
  for( int iEntry=0; iEntry!=nEntries; ++iEntry ) {
    AliAODRecoDecayHF3Prong *dplu = (AliAODRecoDecayHF3Prong*) listDplus->UncheckedAt( iEntry );
    if( !dplu ) continue;
    // APPLYING CUTS
    if( !dplu->HasSelectionBit(AliRDHFCuts::kDplusCuts) ) continue;
    if( !fCutsDStartoKpipi->IsInFiducialAcceptance(dplu->Pt(),dplu->YDplus()) )continue;
    Int_t topCut = fCutsDStartoKpipi->IsSelected( dplu, AliRDHFCuts::kAll );
    if(!topCut) continue;
    Double_t massDp=dplu->InvMassDplus();
    // TO HANDLE AUTOCORRELATIONS
    nIDs[0] = ((AliAODTrack*)dplu->GetDaughter(0))->GetID();
    nIDs[1] = ((AliAODTrack*)dplu->GetDaughter(1))->GetID();
    nIDs[2] = ((AliAODTrack*)dplu->GetDaughter(2))->GetID();
    // ADDING TRACK
    MakeTrack(massDp, dplu->Pt(), dplu->Phi(), dplu->Eta(), 3, nIDs);
    if(fDebugV2) printf("   ᶫInjecting Dplus candidate %d\n",iEntry);
  }
}
//_______________________________________________________________________________
void AliAnalysisTaskFlowD2H::MakeTrack( Double_t mass, Double_t pt, Double_t phi, 
					Double_t eta, Int_t nDau, const Int_t *iID ) {
  // create track for flow tasks
  Bool_t overwrite = kTRUE;
  AliFlowCandidateTrack *oTrack = (static_cast<AliFlowCandidateTrack*> (fCandidates->At( fCandidates->GetLast()+1 )));
  if( !oTrack ) { // creates new
    oTrack = new AliFlowCandidateTrack();
    overwrite = kFALSE;
  } else { // overwrites
    oTrack->ClearMe();
  }
  oTrack->SetMass(mass);
  oTrack->SetPt(pt);
  oTrack->SetPhi(phi);
  oTrack->SetEta(eta);
  for(Int_t iDau=0; iDau!=nDau; ++iDau)
    oTrack->AddDaughter(iID[iDau]);
  oTrack->SetForPOISelection(kTRUE);
  oTrack->SetForRPSelection(kFALSE);
  if(overwrite) {
    fCandidates->SetLast( fCandidates->GetLast()+1 );
  } else {
    fCandidates->AddLast(oTrack);
  }
  return;
}

void AliAnalysisTaskFlowD2H::SetCommonConstants(Int_t massBins, Double_t minMass, Double_t maxMass, Int_t ptWidth) {
  fMassBins = massBins;
  fMinMass = minMass;
  fMaxMass = maxMass;
  fPtBinWidth = ptWidth;
}

//=======================================================================
void AliAnalysisTaskFlowD2H::Terminate(Option_t *)
{
  //terminate
}
 AliAnalysisTaskFlowD2H.cxx:1
 AliAnalysisTaskFlowD2H.cxx:2
 AliAnalysisTaskFlowD2H.cxx:3
 AliAnalysisTaskFlowD2H.cxx:4
 AliAnalysisTaskFlowD2H.cxx:5
 AliAnalysisTaskFlowD2H.cxx:6
 AliAnalysisTaskFlowD2H.cxx:7
 AliAnalysisTaskFlowD2H.cxx:8
 AliAnalysisTaskFlowD2H.cxx:9
 AliAnalysisTaskFlowD2H.cxx:10
 AliAnalysisTaskFlowD2H.cxx:11
 AliAnalysisTaskFlowD2H.cxx:12
 AliAnalysisTaskFlowD2H.cxx:13
 AliAnalysisTaskFlowD2H.cxx:14
 AliAnalysisTaskFlowD2H.cxx:15
 AliAnalysisTaskFlowD2H.cxx:16
 AliAnalysisTaskFlowD2H.cxx:17
 AliAnalysisTaskFlowD2H.cxx:18
 AliAnalysisTaskFlowD2H.cxx:19
 AliAnalysisTaskFlowD2H.cxx:20
 AliAnalysisTaskFlowD2H.cxx:21
 AliAnalysisTaskFlowD2H.cxx:22
 AliAnalysisTaskFlowD2H.cxx:23
 AliAnalysisTaskFlowD2H.cxx:24
 AliAnalysisTaskFlowD2H.cxx:25
 AliAnalysisTaskFlowD2H.cxx:26
 AliAnalysisTaskFlowD2H.cxx:27
 AliAnalysisTaskFlowD2H.cxx:28
 AliAnalysisTaskFlowD2H.cxx:29
 AliAnalysisTaskFlowD2H.cxx:30
 AliAnalysisTaskFlowD2H.cxx:31
 AliAnalysisTaskFlowD2H.cxx:32
 AliAnalysisTaskFlowD2H.cxx:33
 AliAnalysisTaskFlowD2H.cxx:34
 AliAnalysisTaskFlowD2H.cxx:35
 AliAnalysisTaskFlowD2H.cxx:36
 AliAnalysisTaskFlowD2H.cxx:37
 AliAnalysisTaskFlowD2H.cxx:38
 AliAnalysisTaskFlowD2H.cxx:39
 AliAnalysisTaskFlowD2H.cxx:40
 AliAnalysisTaskFlowD2H.cxx:41
 AliAnalysisTaskFlowD2H.cxx:42
 AliAnalysisTaskFlowD2H.cxx:43
 AliAnalysisTaskFlowD2H.cxx:44
 AliAnalysisTaskFlowD2H.cxx:45
 AliAnalysisTaskFlowD2H.cxx:46
 AliAnalysisTaskFlowD2H.cxx:47
 AliAnalysisTaskFlowD2H.cxx:48
 AliAnalysisTaskFlowD2H.cxx:49
 AliAnalysisTaskFlowD2H.cxx:50
 AliAnalysisTaskFlowD2H.cxx:51
 AliAnalysisTaskFlowD2H.cxx:52
 AliAnalysisTaskFlowD2H.cxx:53
 AliAnalysisTaskFlowD2H.cxx:54
 AliAnalysisTaskFlowD2H.cxx:55
 AliAnalysisTaskFlowD2H.cxx:56
 AliAnalysisTaskFlowD2H.cxx:57
 AliAnalysisTaskFlowD2H.cxx:58
 AliAnalysisTaskFlowD2H.cxx:59
 AliAnalysisTaskFlowD2H.cxx:60
 AliAnalysisTaskFlowD2H.cxx:61
 AliAnalysisTaskFlowD2H.cxx:62
 AliAnalysisTaskFlowD2H.cxx:63
 AliAnalysisTaskFlowD2H.cxx:64
 AliAnalysisTaskFlowD2H.cxx:65
 AliAnalysisTaskFlowD2H.cxx:66
 AliAnalysisTaskFlowD2H.cxx:67
 AliAnalysisTaskFlowD2H.cxx:68
 AliAnalysisTaskFlowD2H.cxx:69
 AliAnalysisTaskFlowD2H.cxx:70
 AliAnalysisTaskFlowD2H.cxx:71
 AliAnalysisTaskFlowD2H.cxx:72
 AliAnalysisTaskFlowD2H.cxx:73
 AliAnalysisTaskFlowD2H.cxx:74
 AliAnalysisTaskFlowD2H.cxx:75
 AliAnalysisTaskFlowD2H.cxx:76
 AliAnalysisTaskFlowD2H.cxx:77
 AliAnalysisTaskFlowD2H.cxx:78
 AliAnalysisTaskFlowD2H.cxx:79
 AliAnalysisTaskFlowD2H.cxx:80
 AliAnalysisTaskFlowD2H.cxx:81
 AliAnalysisTaskFlowD2H.cxx:82
 AliAnalysisTaskFlowD2H.cxx:83
 AliAnalysisTaskFlowD2H.cxx:84
 AliAnalysisTaskFlowD2H.cxx:85
 AliAnalysisTaskFlowD2H.cxx:86
 AliAnalysisTaskFlowD2H.cxx:87
 AliAnalysisTaskFlowD2H.cxx:88
 AliAnalysisTaskFlowD2H.cxx:89
 AliAnalysisTaskFlowD2H.cxx:90
 AliAnalysisTaskFlowD2H.cxx:91
 AliAnalysisTaskFlowD2H.cxx:92
 AliAnalysisTaskFlowD2H.cxx:93
 AliAnalysisTaskFlowD2H.cxx:94
 AliAnalysisTaskFlowD2H.cxx:95
 AliAnalysisTaskFlowD2H.cxx:96
 AliAnalysisTaskFlowD2H.cxx:97
 AliAnalysisTaskFlowD2H.cxx:98
 AliAnalysisTaskFlowD2H.cxx:99
 AliAnalysisTaskFlowD2H.cxx:100
 AliAnalysisTaskFlowD2H.cxx:101
 AliAnalysisTaskFlowD2H.cxx:102
 AliAnalysisTaskFlowD2H.cxx:103
 AliAnalysisTaskFlowD2H.cxx:104
 AliAnalysisTaskFlowD2H.cxx:105
 AliAnalysisTaskFlowD2H.cxx:106
 AliAnalysisTaskFlowD2H.cxx:107
 AliAnalysisTaskFlowD2H.cxx:108
 AliAnalysisTaskFlowD2H.cxx:109
 AliAnalysisTaskFlowD2H.cxx:110
 AliAnalysisTaskFlowD2H.cxx:111
 AliAnalysisTaskFlowD2H.cxx:112
 AliAnalysisTaskFlowD2H.cxx:113
 AliAnalysisTaskFlowD2H.cxx:114
 AliAnalysisTaskFlowD2H.cxx:115
 AliAnalysisTaskFlowD2H.cxx:116
 AliAnalysisTaskFlowD2H.cxx:117
 AliAnalysisTaskFlowD2H.cxx:118
 AliAnalysisTaskFlowD2H.cxx:119
 AliAnalysisTaskFlowD2H.cxx:120
 AliAnalysisTaskFlowD2H.cxx:121
 AliAnalysisTaskFlowD2H.cxx:122
 AliAnalysisTaskFlowD2H.cxx:123
 AliAnalysisTaskFlowD2H.cxx:124
 AliAnalysisTaskFlowD2H.cxx:125
 AliAnalysisTaskFlowD2H.cxx:126
 AliAnalysisTaskFlowD2H.cxx:127
 AliAnalysisTaskFlowD2H.cxx:128
 AliAnalysisTaskFlowD2H.cxx:129
 AliAnalysisTaskFlowD2H.cxx:130
 AliAnalysisTaskFlowD2H.cxx:131
 AliAnalysisTaskFlowD2H.cxx:132
 AliAnalysisTaskFlowD2H.cxx:133
 AliAnalysisTaskFlowD2H.cxx:134
 AliAnalysisTaskFlowD2H.cxx:135
 AliAnalysisTaskFlowD2H.cxx:136
 AliAnalysisTaskFlowD2H.cxx:137
 AliAnalysisTaskFlowD2H.cxx:138
 AliAnalysisTaskFlowD2H.cxx:139
 AliAnalysisTaskFlowD2H.cxx:140
 AliAnalysisTaskFlowD2H.cxx:141
 AliAnalysisTaskFlowD2H.cxx:142
 AliAnalysisTaskFlowD2H.cxx:143
 AliAnalysisTaskFlowD2H.cxx:144
 AliAnalysisTaskFlowD2H.cxx:145
 AliAnalysisTaskFlowD2H.cxx:146
 AliAnalysisTaskFlowD2H.cxx:147
 AliAnalysisTaskFlowD2H.cxx:148
 AliAnalysisTaskFlowD2H.cxx:149
 AliAnalysisTaskFlowD2H.cxx:150
 AliAnalysisTaskFlowD2H.cxx:151
 AliAnalysisTaskFlowD2H.cxx:152
 AliAnalysisTaskFlowD2H.cxx:153
 AliAnalysisTaskFlowD2H.cxx:154
 AliAnalysisTaskFlowD2H.cxx:155
 AliAnalysisTaskFlowD2H.cxx:156
 AliAnalysisTaskFlowD2H.cxx:157
 AliAnalysisTaskFlowD2H.cxx:158
 AliAnalysisTaskFlowD2H.cxx:159
 AliAnalysisTaskFlowD2H.cxx:160
 AliAnalysisTaskFlowD2H.cxx:161
 AliAnalysisTaskFlowD2H.cxx:162
 AliAnalysisTaskFlowD2H.cxx:163
 AliAnalysisTaskFlowD2H.cxx:164
 AliAnalysisTaskFlowD2H.cxx:165
 AliAnalysisTaskFlowD2H.cxx:166
 AliAnalysisTaskFlowD2H.cxx:167
 AliAnalysisTaskFlowD2H.cxx:168
 AliAnalysisTaskFlowD2H.cxx:169
 AliAnalysisTaskFlowD2H.cxx:170
 AliAnalysisTaskFlowD2H.cxx:171
 AliAnalysisTaskFlowD2H.cxx:172
 AliAnalysisTaskFlowD2H.cxx:173
 AliAnalysisTaskFlowD2H.cxx:174
 AliAnalysisTaskFlowD2H.cxx:175
 AliAnalysisTaskFlowD2H.cxx:176
 AliAnalysisTaskFlowD2H.cxx:177
 AliAnalysisTaskFlowD2H.cxx:178
 AliAnalysisTaskFlowD2H.cxx:179
 AliAnalysisTaskFlowD2H.cxx:180
 AliAnalysisTaskFlowD2H.cxx:181
 AliAnalysisTaskFlowD2H.cxx:182
 AliAnalysisTaskFlowD2H.cxx:183
 AliAnalysisTaskFlowD2H.cxx:184
 AliAnalysisTaskFlowD2H.cxx:185
 AliAnalysisTaskFlowD2H.cxx:186
 AliAnalysisTaskFlowD2H.cxx:187
 AliAnalysisTaskFlowD2H.cxx:188
 AliAnalysisTaskFlowD2H.cxx:189
 AliAnalysisTaskFlowD2H.cxx:190
 AliAnalysisTaskFlowD2H.cxx:191
 AliAnalysisTaskFlowD2H.cxx:192
 AliAnalysisTaskFlowD2H.cxx:193
 AliAnalysisTaskFlowD2H.cxx:194
 AliAnalysisTaskFlowD2H.cxx:195
 AliAnalysisTaskFlowD2H.cxx:196
 AliAnalysisTaskFlowD2H.cxx:197
 AliAnalysisTaskFlowD2H.cxx:198
 AliAnalysisTaskFlowD2H.cxx:199
 AliAnalysisTaskFlowD2H.cxx:200
 AliAnalysisTaskFlowD2H.cxx:201
 AliAnalysisTaskFlowD2H.cxx:202
 AliAnalysisTaskFlowD2H.cxx:203
 AliAnalysisTaskFlowD2H.cxx:204
 AliAnalysisTaskFlowD2H.cxx:205
 AliAnalysisTaskFlowD2H.cxx:206
 AliAnalysisTaskFlowD2H.cxx:207
 AliAnalysisTaskFlowD2H.cxx:208
 AliAnalysisTaskFlowD2H.cxx:209
 AliAnalysisTaskFlowD2H.cxx:210
 AliAnalysisTaskFlowD2H.cxx:211
 AliAnalysisTaskFlowD2H.cxx:212
 AliAnalysisTaskFlowD2H.cxx:213
 AliAnalysisTaskFlowD2H.cxx:214
 AliAnalysisTaskFlowD2H.cxx:215
 AliAnalysisTaskFlowD2H.cxx:216
 AliAnalysisTaskFlowD2H.cxx:217
 AliAnalysisTaskFlowD2H.cxx:218
 AliAnalysisTaskFlowD2H.cxx:219
 AliAnalysisTaskFlowD2H.cxx:220
 AliAnalysisTaskFlowD2H.cxx:221
 AliAnalysisTaskFlowD2H.cxx:222
 AliAnalysisTaskFlowD2H.cxx:223
 AliAnalysisTaskFlowD2H.cxx:224
 AliAnalysisTaskFlowD2H.cxx:225
 AliAnalysisTaskFlowD2H.cxx:226
 AliAnalysisTaskFlowD2H.cxx:227
 AliAnalysisTaskFlowD2H.cxx:228
 AliAnalysisTaskFlowD2H.cxx:229
 AliAnalysisTaskFlowD2H.cxx:230
 AliAnalysisTaskFlowD2H.cxx:231
 AliAnalysisTaskFlowD2H.cxx:232
 AliAnalysisTaskFlowD2H.cxx:233
 AliAnalysisTaskFlowD2H.cxx:234
 AliAnalysisTaskFlowD2H.cxx:235
 AliAnalysisTaskFlowD2H.cxx:236
 AliAnalysisTaskFlowD2H.cxx:237
 AliAnalysisTaskFlowD2H.cxx:238
 AliAnalysisTaskFlowD2H.cxx:239
 AliAnalysisTaskFlowD2H.cxx:240
 AliAnalysisTaskFlowD2H.cxx:241
 AliAnalysisTaskFlowD2H.cxx:242
 AliAnalysisTaskFlowD2H.cxx:243
 AliAnalysisTaskFlowD2H.cxx:244
 AliAnalysisTaskFlowD2H.cxx:245
 AliAnalysisTaskFlowD2H.cxx:246
 AliAnalysisTaskFlowD2H.cxx:247
 AliAnalysisTaskFlowD2H.cxx:248
 AliAnalysisTaskFlowD2H.cxx:249
 AliAnalysisTaskFlowD2H.cxx:250
 AliAnalysisTaskFlowD2H.cxx:251
 AliAnalysisTaskFlowD2H.cxx:252
 AliAnalysisTaskFlowD2H.cxx:253
 AliAnalysisTaskFlowD2H.cxx:254
 AliAnalysisTaskFlowD2H.cxx:255
 AliAnalysisTaskFlowD2H.cxx:256
 AliAnalysisTaskFlowD2H.cxx:257
 AliAnalysisTaskFlowD2H.cxx:258
 AliAnalysisTaskFlowD2H.cxx:259
 AliAnalysisTaskFlowD2H.cxx:260
 AliAnalysisTaskFlowD2H.cxx:261
 AliAnalysisTaskFlowD2H.cxx:262
 AliAnalysisTaskFlowD2H.cxx:263
 AliAnalysisTaskFlowD2H.cxx:264
 AliAnalysisTaskFlowD2H.cxx:265
 AliAnalysisTaskFlowD2H.cxx:266
 AliAnalysisTaskFlowD2H.cxx:267
 AliAnalysisTaskFlowD2H.cxx:268
 AliAnalysisTaskFlowD2H.cxx:269
 AliAnalysisTaskFlowD2H.cxx:270
 AliAnalysisTaskFlowD2H.cxx:271
 AliAnalysisTaskFlowD2H.cxx:272
 AliAnalysisTaskFlowD2H.cxx:273
 AliAnalysisTaskFlowD2H.cxx:274
 AliAnalysisTaskFlowD2H.cxx:275
 AliAnalysisTaskFlowD2H.cxx:276
 AliAnalysisTaskFlowD2H.cxx:277
 AliAnalysisTaskFlowD2H.cxx:278
 AliAnalysisTaskFlowD2H.cxx:279
 AliAnalysisTaskFlowD2H.cxx:280
 AliAnalysisTaskFlowD2H.cxx:281
 AliAnalysisTaskFlowD2H.cxx:282
 AliAnalysisTaskFlowD2H.cxx:283
 AliAnalysisTaskFlowD2H.cxx:284
 AliAnalysisTaskFlowD2H.cxx:285
 AliAnalysisTaskFlowD2H.cxx:286
 AliAnalysisTaskFlowD2H.cxx:287
 AliAnalysisTaskFlowD2H.cxx:288
 AliAnalysisTaskFlowD2H.cxx:289
 AliAnalysisTaskFlowD2H.cxx:290
 AliAnalysisTaskFlowD2H.cxx:291
 AliAnalysisTaskFlowD2H.cxx:292
 AliAnalysisTaskFlowD2H.cxx:293
 AliAnalysisTaskFlowD2H.cxx:294
 AliAnalysisTaskFlowD2H.cxx:295
 AliAnalysisTaskFlowD2H.cxx:296
 AliAnalysisTaskFlowD2H.cxx:297
 AliAnalysisTaskFlowD2H.cxx:298
 AliAnalysisTaskFlowD2H.cxx:299
 AliAnalysisTaskFlowD2H.cxx:300
 AliAnalysisTaskFlowD2H.cxx:301
 AliAnalysisTaskFlowD2H.cxx:302
 AliAnalysisTaskFlowD2H.cxx:303
 AliAnalysisTaskFlowD2H.cxx:304
 AliAnalysisTaskFlowD2H.cxx:305
 AliAnalysisTaskFlowD2H.cxx:306
 AliAnalysisTaskFlowD2H.cxx:307
 AliAnalysisTaskFlowD2H.cxx:308
 AliAnalysisTaskFlowD2H.cxx:309
 AliAnalysisTaskFlowD2H.cxx:310
 AliAnalysisTaskFlowD2H.cxx:311
 AliAnalysisTaskFlowD2H.cxx:312
 AliAnalysisTaskFlowD2H.cxx:313
 AliAnalysisTaskFlowD2H.cxx:314
 AliAnalysisTaskFlowD2H.cxx:315
 AliAnalysisTaskFlowD2H.cxx:316
 AliAnalysisTaskFlowD2H.cxx:317
 AliAnalysisTaskFlowD2H.cxx:318
 AliAnalysisTaskFlowD2H.cxx:319
 AliAnalysisTaskFlowD2H.cxx:320
 AliAnalysisTaskFlowD2H.cxx:321
 AliAnalysisTaskFlowD2H.cxx:322
 AliAnalysisTaskFlowD2H.cxx:323
 AliAnalysisTaskFlowD2H.cxx:324
 AliAnalysisTaskFlowD2H.cxx:325
 AliAnalysisTaskFlowD2H.cxx:326
 AliAnalysisTaskFlowD2H.cxx:327
 AliAnalysisTaskFlowD2H.cxx:328
 AliAnalysisTaskFlowD2H.cxx:329
 AliAnalysisTaskFlowD2H.cxx:330
 AliAnalysisTaskFlowD2H.cxx:331
 AliAnalysisTaskFlowD2H.cxx:332
 AliAnalysisTaskFlowD2H.cxx:333
 AliAnalysisTaskFlowD2H.cxx:334
 AliAnalysisTaskFlowD2H.cxx:335
 AliAnalysisTaskFlowD2H.cxx:336
 AliAnalysisTaskFlowD2H.cxx:337
 AliAnalysisTaskFlowD2H.cxx:338
 AliAnalysisTaskFlowD2H.cxx:339
 AliAnalysisTaskFlowD2H.cxx:340
 AliAnalysisTaskFlowD2H.cxx:341
 AliAnalysisTaskFlowD2H.cxx:342
 AliAnalysisTaskFlowD2H.cxx:343
 AliAnalysisTaskFlowD2H.cxx:344
 AliAnalysisTaskFlowD2H.cxx:345
 AliAnalysisTaskFlowD2H.cxx:346
 AliAnalysisTaskFlowD2H.cxx:347
 AliAnalysisTaskFlowD2H.cxx:348
 AliAnalysisTaskFlowD2H.cxx:349
 AliAnalysisTaskFlowD2H.cxx:350
 AliAnalysisTaskFlowD2H.cxx:351
 AliAnalysisTaskFlowD2H.cxx:352
 AliAnalysisTaskFlowD2H.cxx:353
 AliAnalysisTaskFlowD2H.cxx:354
 AliAnalysisTaskFlowD2H.cxx:355
 AliAnalysisTaskFlowD2H.cxx:356
 AliAnalysisTaskFlowD2H.cxx:357
 AliAnalysisTaskFlowD2H.cxx:358
 AliAnalysisTaskFlowD2H.cxx:359
 AliAnalysisTaskFlowD2H.cxx:360
 AliAnalysisTaskFlowD2H.cxx:361
 AliAnalysisTaskFlowD2H.cxx:362
 AliAnalysisTaskFlowD2H.cxx:363
 AliAnalysisTaskFlowD2H.cxx:364
 AliAnalysisTaskFlowD2H.cxx:365
 AliAnalysisTaskFlowD2H.cxx:366
 AliAnalysisTaskFlowD2H.cxx:367
 AliAnalysisTaskFlowD2H.cxx:368
 AliAnalysisTaskFlowD2H.cxx:369
 AliAnalysisTaskFlowD2H.cxx:370
 AliAnalysisTaskFlowD2H.cxx:371
 AliAnalysisTaskFlowD2H.cxx:372
 AliAnalysisTaskFlowD2H.cxx:373
 AliAnalysisTaskFlowD2H.cxx:374
 AliAnalysisTaskFlowD2H.cxx:375
 AliAnalysisTaskFlowD2H.cxx:376
 AliAnalysisTaskFlowD2H.cxx:377
 AliAnalysisTaskFlowD2H.cxx:378
 AliAnalysisTaskFlowD2H.cxx:379
 AliAnalysisTaskFlowD2H.cxx:380
 AliAnalysisTaskFlowD2H.cxx:381
 AliAnalysisTaskFlowD2H.cxx:382