ROOT logo
//DEFINITION OF A FEW CONSTANTS
const Double_t ptmin  =   0.0 ;
const Double_t ptmax  =   1.0 ;
const Double_t etamin =  -1.5 ;
const Double_t etamax =   1.5 ;
const Int_t    PDG    =    211; 
const Int_t    minclustersTPC = 40 ;
//----------------------------------------------------

/*

Macro to prepare the necessary objects to perform spectrum unfolding 
(to be used by the class AliCFUnfolding)

Note that the efficiency calculation (and therfore the container filling) 
has to be done using MC values (and not reconstructed values)

This macro creates the following objects :
- (AliCFContainer) container    : used to calculate the efficiency (see AliCFEffGrid)
- (THnSparseD)     correlation  : this is the response matrix

Once you have run this macro, you may launch unfolding using 
the example macro CORRFW/test/testUnfolding.C

*/

Bool_t AliCFTaskForUnfolding()
{
  
  TBenchmark benchmark;
  benchmark.Start("AliSingleTrackTask");

  AliLog::SetGlobalDebugLevel(0);

  Load() ; //load the required libraries

  TChain * analysisChain ;
  printf("\n\nRunning on local file, please check the path\n\n");

  analysisChain = new TChain("esdTree");
  analysisChain->Add("your_data_path/001/AliESDs.root");
  analysisChain->Add("your_data_path/002/AliESDs.root");

  
  Info("AliCFTaskForUnfolding",Form("CHAIN HAS %d ENTRIES",(Int_t)analysisChain->GetEntries()));

  //CONTAINER DEFINITION
  Info("AliCFTaskForUnfolding","SETUP CONTAINER");
  //the sensitive variables (2 in this example), their indices
  UInt_t ipt = 0;
  UInt_t ieta  = 1;
  //Setting up the container grid... 
  UInt_t nstep = 3 ; //number of selection steps MC 
  const Int_t nvar   = 2 ; //number of variables on the grid:pt,eta
  const Int_t nbin[nvar] = {20,20} ;

  //arrays for the number of bins in each dimension
  Int_t iBin[nvar];
  iBin[0]=nbin[0];
  iBin[1]=nbin[1];

  //arrays for lower bounds :
  Double_t *binLim0=new Double_t[nbin[0]+1];
  Double_t *binLim1=new Double_t[nbin[1]+1];

  //values for bin lower bounds
  for(Int_t i=0; i<=nbin[0]; i++) binLim0[i]=(Double_t)ptmin  + ( ptmax- ptmin)/nbin[0]*(Double_t)i ; 
  for(Int_t i=0; i<=nbin[1]; i++) binLim1[i]=(Double_t)etamin + (etamax-etamin)/nbin[1]*(Double_t)i ; 

  //one "container" for MC
  AliCFContainer* container = new AliCFContainer("container","container for tracks",nstep,nvar,iBin);
  //setting the bin limits
  container -> SetBinLimits(ipt,binLim0);
  container -> SetBinLimits(ieta ,binLim1);
  container -> SetVarTitle(ipt,"pT");
  container -> SetVarTitle(ieta,"#eta");

  // Gen-Level kinematic cuts
  AliCFTrackKineCuts *mcKineCuts = new AliCFTrackKineCuts("mcKineCuts","MC-level kinematic cuts");
  mcKineCuts->SetPtRange (ptmin ,ptmax );
  mcKineCuts->SetEtaRange(etamin,etamax);

  //Particle-Level cuts:  
  AliCFParticleGenCuts* mcGenCuts = new AliCFParticleGenCuts("mcGenCuts","MC particle generation cuts");
  mcGenCuts->SetRequireIsPrimary();
  mcGenCuts->SetRequirePdgCode(PDG);

  // Rec-Level kinematic cuts
  AliCFTrackKineCuts *recKineCuts = new AliCFTrackKineCuts("recKineCuts","rec-level kine cuts");
  recKineCuts->SetPtRange( ptmin, ptmax);
  recKineCuts->SetPtRange(etamin,etamax);

  AliCFTrackQualityCuts *recQualityCuts = new AliCFTrackQualityCuts("recQualityCuts","rec-level quality cuts");
  recQualityCuts->SetMinNClusterTPC(minclustersTPC);

  AliCFTrackIsPrimaryCuts *recIsPrimaryCuts = new AliCFTrackIsPrimaryCuts("recIsPrimaryCuts","rec-level isPrimary cuts");
  recIsPrimaryCuts->SetMaxNSigmaToVertex(3);

  printf("CREATE MC KINE CUTS\n");
  TObjArray* mcList = new TObjArray(0) ;
  mcList->AddLast(mcKineCuts);
  mcList->AddLast(mcGenCuts);

  printf("CREATE RECONSTRUCTION CUTS\n");
  TObjArray* recList = new TObjArray(0) ;
  recList->AddLast(recKineCuts);
  recList->AddLast(recQualityCuts);
  recList->AddLast(recIsPrimaryCuts);

  TObjArray* emptyList = new TObjArray(0);

  //CREATE THE INTERFACE TO CORRECTION FRAMEWORK USED IN THE TASK
  printf("CREATE INTERFACE AND CUTS\n");
  AliCFManager* man = new AliCFManager() ;
  man->SetNStepEvent(0);
  man->SetParticleContainer(container);
  man->SetParticleCutsList(0,mcList);
  man->SetParticleCutsList(1,recList);
  man->SetParticleCutsList(2,emptyList); // this is special step for monte carlo spectrum

  //CREATE THE TASK
  printf("CREATE TASK\n");
  // create the task
  AliCFTaskForUnfolding *task = new AliCFTaskForUnfolding("AliCFTaskForUnfolding");
  task->SetCFManager(man); //here is set the CF manager

  //create correlation matrix for unfolding
  Int_t thnDim[2*nvar];
  for (int k=0; k<nvar; k++) {
    //first half  : reconstructed 
    //second half : MC
    thnDim[k]      = nbin[k];
    thnDim[k+nvar] = nbin[k];
  }
  THnSparseD* correlation = new THnSparseD("correlation","THnSparse with correlations",2*nvar,thnDim);
  Double_t** binEdges = new Double_t[nvar];
  for (int k=0; k<nvar; k++) {
    binEdges[k]=new Double_t[nbin[k]+1];
    container->GetBinLimits(k,binEdges[k]);
    correlation->SetBinEdges(k,binEdges[k]);
    correlation->SetBinEdges(k+nvar,binEdges[k]);
  }
  correlation->Sumw2();
  task->SetCorrelationMatrix(correlation);
  //---

  //SETUP THE ANALYSIS MANAGER TO READ INPUT CHAIN AND WRITE DESIRED OUTPUTS
  printf("CREATE ANALYSIS MANAGER\n");
  // Make the analysis manager
  AliAnalysisManager *mgr = new AliAnalysisManager("TestManager");
  mgr->SetAnalysisType(AliAnalysisManager::kLocalAnalysis);

  AliMCEventHandler*  mcHandler = new AliMCEventHandler();
  mgr->SetMCtruthEventHandler(mcHandler);
 
  AliInputEventHandler* dataHandler = new AliESDInputHandler();
  mgr->SetInputEventHandler(dataHandler);

  // Create and connect containers for input/output

  //------ input data ------
  AliAnalysisDataContainer *cinput0  = mgr->CreateContainer("cchain0",TChain::Class(),AliAnalysisManager::kInputContainer);

  // ----- output data -----
  
  //slot 0 : default output tree (by default handled by AliAnalysisTaskSE)
  AliAnalysisDataContainer *coutput0 = mgr->CreateContainer("ctree0", TTree::Class(),AliAnalysisManager::kOutputContainer,"output.root");

  //now comes user's output objects :
  
  // output TH1I for event counting
  AliAnalysisDataContainer *coutput1 = mgr->CreateContainer("chist0", TH1I::Class(),AliAnalysisManager::kOutputContainer,"output.root");
  // output Correction Framework Container (for acceptance & efficiency calculations)
  AliAnalysisDataContainer *coutput2 = mgr->CreateContainer("ccontainer0", AliCFContainer::Class(),AliAnalysisManager::kOutputContainer,"output.root");
  AliAnalysisDataContainer *coutput3 = mgr->CreateContainer("corr0", THnSparseD::Class(),AliAnalysisManager::kOutputContainer,"output.root");

  cinput0->SetData(analysisChain);

  mgr->AddTask(task);
  mgr->ConnectInput(task,0,mgr->GetCommonInputContainer());
  mgr->ConnectOutput(task,0,coutput0);
  mgr->ConnectOutput(task,1,coutput1);
  mgr->ConnectOutput(task,2,coutput2);
  mgr->ConnectOutput(task,3,coutput3);
 
  printf("READY TO RUN\n");
  //RUN !!!
  if (mgr->InitAnalysis()) {
    mgr->PrintStatus();
    mgr->StartAnalysis("local",analysisChain);
  }

  benchmark.Stop("AliSingleTrackTask");
  benchmark.Show("AliSingleTrackTask");

  return kTRUE ;
}

void Load() {

  //load the required aliroot libraries
  gSystem->Load("libANALYSIS") ;
  gSystem->Load("libANALYSISalice") ;
  gSystem->Load("libCORRFW.so") ;

  //compile online the task class
  gSystem->SetIncludePath("-I. -I$ALICE_ROOT/include -I$ROOTSYS/include");
  gROOT->LoadMacro("./AliCFTaskForUnfolding.cxx+");
}
 AliCFTaskForUnfolding.C:1
 AliCFTaskForUnfolding.C:2
 AliCFTaskForUnfolding.C:3
 AliCFTaskForUnfolding.C:4
 AliCFTaskForUnfolding.C:5
 AliCFTaskForUnfolding.C:6
 AliCFTaskForUnfolding.C:7
 AliCFTaskForUnfolding.C:8
 AliCFTaskForUnfolding.C:9
 AliCFTaskForUnfolding.C:10
 AliCFTaskForUnfolding.C:11
 AliCFTaskForUnfolding.C:12
 AliCFTaskForUnfolding.C:13
 AliCFTaskForUnfolding.C:14
 AliCFTaskForUnfolding.C:15
 AliCFTaskForUnfolding.C:16
 AliCFTaskForUnfolding.C:17
 AliCFTaskForUnfolding.C:18
 AliCFTaskForUnfolding.C:19
 AliCFTaskForUnfolding.C:20
 AliCFTaskForUnfolding.C:21
 AliCFTaskForUnfolding.C:22
 AliCFTaskForUnfolding.C:23
 AliCFTaskForUnfolding.C:24
 AliCFTaskForUnfolding.C:25
 AliCFTaskForUnfolding.C:26
 AliCFTaskForUnfolding.C:27
 AliCFTaskForUnfolding.C:28
 AliCFTaskForUnfolding.C:29
 AliCFTaskForUnfolding.C:30
 AliCFTaskForUnfolding.C:31
 AliCFTaskForUnfolding.C:32
 AliCFTaskForUnfolding.C:33
 AliCFTaskForUnfolding.C:34
 AliCFTaskForUnfolding.C:35
 AliCFTaskForUnfolding.C:36
 AliCFTaskForUnfolding.C:37
 AliCFTaskForUnfolding.C:38
 AliCFTaskForUnfolding.C:39
 AliCFTaskForUnfolding.C:40
 AliCFTaskForUnfolding.C:41
 AliCFTaskForUnfolding.C:42
 AliCFTaskForUnfolding.C:43
 AliCFTaskForUnfolding.C:44
 AliCFTaskForUnfolding.C:45
 AliCFTaskForUnfolding.C:46
 AliCFTaskForUnfolding.C:47
 AliCFTaskForUnfolding.C:48
 AliCFTaskForUnfolding.C:49
 AliCFTaskForUnfolding.C:50
 AliCFTaskForUnfolding.C:51
 AliCFTaskForUnfolding.C:52
 AliCFTaskForUnfolding.C:53
 AliCFTaskForUnfolding.C:54
 AliCFTaskForUnfolding.C:55
 AliCFTaskForUnfolding.C:56
 AliCFTaskForUnfolding.C:57
 AliCFTaskForUnfolding.C:58
 AliCFTaskForUnfolding.C:59
 AliCFTaskForUnfolding.C:60
 AliCFTaskForUnfolding.C:61
 AliCFTaskForUnfolding.C:62
 AliCFTaskForUnfolding.C:63
 AliCFTaskForUnfolding.C:64
 AliCFTaskForUnfolding.C:65
 AliCFTaskForUnfolding.C:66
 AliCFTaskForUnfolding.C:67
 AliCFTaskForUnfolding.C:68
 AliCFTaskForUnfolding.C:69
 AliCFTaskForUnfolding.C:70
 AliCFTaskForUnfolding.C:71
 AliCFTaskForUnfolding.C:72
 AliCFTaskForUnfolding.C:73
 AliCFTaskForUnfolding.C:74
 AliCFTaskForUnfolding.C:75
 AliCFTaskForUnfolding.C:76
 AliCFTaskForUnfolding.C:77
 AliCFTaskForUnfolding.C:78
 AliCFTaskForUnfolding.C:79
 AliCFTaskForUnfolding.C:80
 AliCFTaskForUnfolding.C:81
 AliCFTaskForUnfolding.C:82
 AliCFTaskForUnfolding.C:83
 AliCFTaskForUnfolding.C:84
 AliCFTaskForUnfolding.C:85
 AliCFTaskForUnfolding.C:86
 AliCFTaskForUnfolding.C:87
 AliCFTaskForUnfolding.C:88
 AliCFTaskForUnfolding.C:89
 AliCFTaskForUnfolding.C:90
 AliCFTaskForUnfolding.C:91
 AliCFTaskForUnfolding.C:92
 AliCFTaskForUnfolding.C:93
 AliCFTaskForUnfolding.C:94
 AliCFTaskForUnfolding.C:95
 AliCFTaskForUnfolding.C:96
 AliCFTaskForUnfolding.C:97
 AliCFTaskForUnfolding.C:98
 AliCFTaskForUnfolding.C:99
 AliCFTaskForUnfolding.C:100
 AliCFTaskForUnfolding.C:101
 AliCFTaskForUnfolding.C:102
 AliCFTaskForUnfolding.C:103
 AliCFTaskForUnfolding.C:104
 AliCFTaskForUnfolding.C:105
 AliCFTaskForUnfolding.C:106
 AliCFTaskForUnfolding.C:107
 AliCFTaskForUnfolding.C:108
 AliCFTaskForUnfolding.C:109
 AliCFTaskForUnfolding.C:110
 AliCFTaskForUnfolding.C:111
 AliCFTaskForUnfolding.C:112
 AliCFTaskForUnfolding.C:113
 AliCFTaskForUnfolding.C:114
 AliCFTaskForUnfolding.C:115
 AliCFTaskForUnfolding.C:116
 AliCFTaskForUnfolding.C:117
 AliCFTaskForUnfolding.C:118
 AliCFTaskForUnfolding.C:119
 AliCFTaskForUnfolding.C:120
 AliCFTaskForUnfolding.C:121
 AliCFTaskForUnfolding.C:122
 AliCFTaskForUnfolding.C:123
 AliCFTaskForUnfolding.C:124
 AliCFTaskForUnfolding.C:125
 AliCFTaskForUnfolding.C:126
 AliCFTaskForUnfolding.C:127
 AliCFTaskForUnfolding.C:128
 AliCFTaskForUnfolding.C:129
 AliCFTaskForUnfolding.C:130
 AliCFTaskForUnfolding.C:131
 AliCFTaskForUnfolding.C:132
 AliCFTaskForUnfolding.C:133
 AliCFTaskForUnfolding.C:134
 AliCFTaskForUnfolding.C:135
 AliCFTaskForUnfolding.C:136
 AliCFTaskForUnfolding.C:137
 AliCFTaskForUnfolding.C:138
 AliCFTaskForUnfolding.C:139
 AliCFTaskForUnfolding.C:140
 AliCFTaskForUnfolding.C:141
 AliCFTaskForUnfolding.C:142
 AliCFTaskForUnfolding.C:143
 AliCFTaskForUnfolding.C:144
 AliCFTaskForUnfolding.C:145
 AliCFTaskForUnfolding.C:146
 AliCFTaskForUnfolding.C:147
 AliCFTaskForUnfolding.C:148
 AliCFTaskForUnfolding.C:149
 AliCFTaskForUnfolding.C:150
 AliCFTaskForUnfolding.C:151
 AliCFTaskForUnfolding.C:152
 AliCFTaskForUnfolding.C:153
 AliCFTaskForUnfolding.C:154
 AliCFTaskForUnfolding.C:155
 AliCFTaskForUnfolding.C:156
 AliCFTaskForUnfolding.C:157
 AliCFTaskForUnfolding.C:158
 AliCFTaskForUnfolding.C:159
 AliCFTaskForUnfolding.C:160
 AliCFTaskForUnfolding.C:161
 AliCFTaskForUnfolding.C:162
 AliCFTaskForUnfolding.C:163
 AliCFTaskForUnfolding.C:164
 AliCFTaskForUnfolding.C:165
 AliCFTaskForUnfolding.C:166
 AliCFTaskForUnfolding.C:167
 AliCFTaskForUnfolding.C:168
 AliCFTaskForUnfolding.C:169
 AliCFTaskForUnfolding.C:170
 AliCFTaskForUnfolding.C:171
 AliCFTaskForUnfolding.C:172
 AliCFTaskForUnfolding.C:173
 AliCFTaskForUnfolding.C:174
 AliCFTaskForUnfolding.C:175
 AliCFTaskForUnfolding.C:176
 AliCFTaskForUnfolding.C:177
 AliCFTaskForUnfolding.C:178
 AliCFTaskForUnfolding.C:179
 AliCFTaskForUnfolding.C:180
 AliCFTaskForUnfolding.C:181
 AliCFTaskForUnfolding.C:182
 AliCFTaskForUnfolding.C:183
 AliCFTaskForUnfolding.C:184
 AliCFTaskForUnfolding.C:185
 AliCFTaskForUnfolding.C:186
 AliCFTaskForUnfolding.C:187
 AliCFTaskForUnfolding.C:188
 AliCFTaskForUnfolding.C:189
 AliCFTaskForUnfolding.C:190
 AliCFTaskForUnfolding.C:191
 AliCFTaskForUnfolding.C:192
 AliCFTaskForUnfolding.C:193
 AliCFTaskForUnfolding.C:194
 AliCFTaskForUnfolding.C:195
 AliCFTaskForUnfolding.C:196
 AliCFTaskForUnfolding.C:197
 AliCFTaskForUnfolding.C:198
 AliCFTaskForUnfolding.C:199
 AliCFTaskForUnfolding.C:200
 AliCFTaskForUnfolding.C:201
 AliCFTaskForUnfolding.C:202
 AliCFTaskForUnfolding.C:203
 AliCFTaskForUnfolding.C:204
 AliCFTaskForUnfolding.C:205
 AliCFTaskForUnfolding.C:206
 AliCFTaskForUnfolding.C:207
 AliCFTaskForUnfolding.C:208
 AliCFTaskForUnfolding.C:209
 AliCFTaskForUnfolding.C:210