]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/macros/AddTaskDijetHadron.C
Add new DijetHadron task from Taiyo
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / macros / AddTaskDijetHadron.C
1 AliAnalysisTaskDijetHadron* AddTaskDijetHadron(
2   const char *ntracks            = "Tracks",
3   const char *nclusters          = "CaloClusters",
4   const char *njets              = "Jets",
5   const char *nMCtracks            = "TracksMC",
6   const char *nMCclusters          = "CaloClustersMC",
7   const char *nMCjets              = "JetsMC",
8   const char *nembtracks         = "TracksEmbedded",
9   const char *nembclusters       = "CaloClustersEmbedded",
10   const char *nembjets           = "EmbJets",
11   const char *nrandtracks        = "TracksRandomized",
12   const char *nrandclusters      = "CaloClustersRandomized",
13   const char *nPbPbrho               = "Rho",
14   const char *nMCrho               = "RhoMC",
15   const char *nEMBrho               = "RhoEMB",
16   Double_t    jetradius          = 0.2,
17   Double_t    jetareacut         = 0.557,
18   Double_t    trackptcut         = 0.15,
19   Double_t    clusptcut          = 0.30,
20   const char *type               = "TPC",
21   const char *taskname           = "AliAnalysisTaskDijetHadron"
22 )
23 {  
24   // Get the pointer to the existing analysis manager via the static access method.
25   //==============================================================================
26   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
27   if (!mgr)
28   {
29     ::Error("AddTaskHJetCorr", "No analysis manager to connect to.");
30     return NULL;
31   }  
32   
33   // Check the analysis type using the event handlers connected to the analysis manager.
34   //==============================================================================
35   if (!mgr->GetInputEventHandler())
36   {
37     ::Error("AddTaskHJetCorr", "This task requires an input event handler");
38     return NULL;
39   }
40   
41   //-------------------------------------------------------
42   // Init the task and do settings
43   //-------------------------------------------------------
44   TString name;
45   if (strcmp(ntracks, "") == 0 && strcmp(nclusters, "") == 0) 
46     name = Form("%s_%s_R0%d_%s",taskname,nPbPbrho,(Int_t)floor(jetradius*100+0.5),type);
47   else if (strcmp(ntracks, "") == 0) 
48     name = Form("%s_%s_%s_R0%d_%s",taskname,nclusters,nPbPbrho,(Int_t)floor(jetradius*100+0.5),type);
49   else if (strcmp(nclusters, "") == 0) 
50     name = Form("%s_%s_%s_R0%d_%s",taskname,ntracks,nPbPbrho,(Int_t)floor(jetradius*100+0.5),type);
51   else
52     name = Form("%s_%s_%s_%s_R0%d_%s",taskname,ntracks,nclusters,nPbPbrho,(Int_t)floor(jetradius*100+0.5),type);
53
54   AliAnalysisTaskDijetHadron* jetTask = new AliAnalysisTaskDijetHadron(name);
55   jetTask->SetConeRadius(jetradius);
56   jetTask->SetRhoName(nPbPbrho,-1);
57   if (strcmp(type,"TPC")==0) 
58     jetTask->SetConeEtaPhiTPC();
59   else if (strcmp(type,"EMCAL")==0) 
60     jetTask->SetConeEtaPhiEMCAL();
61
62   AliParticleContainer *partCont = jetTask->AddParticleContainer(ntracks);
63   if (partCont) {
64     partCont->SetName("Tracks");
65     partCont->SetParticlePtCut(trackptcut);
66   }
67
68   AliClusterContainer *clusCont = jetTask->AddClusterContainer(nclusters);
69   if (clusCont) {
70     clusCont->SetName("CaloClusters");
71     clusCont->SetClusPtCut(clusptcut);
72   }
73
74   AliJetContainer *jetCont = jetTask->AddJetContainer(njets,type,jetradius);
75   if (jetCont) {
76     jetCont->SetName("Jets");
77     jetCont->SetPercAreaCut(jetareacut);
78     jetCont->SetRhoName(nPbPbrho);
79     jetCont->ConnectParticleContainer(partCont);
80     jetCont->ConnectClusterContainer(clusCont);
81   }
82
83   AliParticleContainer *MCpartCont = jetTask->AddParticleContainer(nMCtracks);
84   if (partCont) {
85     MCpartCont->SetName("MCTracks");
86     MCpartCont->SetParticlePtCut(trackptcut);
87   }
88
89   AliClusterContainer *MCclusCont = jetTask->AddClusterContainer(nMCclusters);
90   if (clusCont) {
91     MCclusCont->SetName("MCCaloClusters");
92     MCclusCont->SetClusPtCut(clusptcut);
93   }
94
95   AliJetContainer *MCjetCont = jetTask->AddJetContainer(nMCjets,type,jetradius);
96   if (jetCont) {
97     MCjetCont->SetName("MCJets");
98     MCjetCont->SetPercAreaCut(jetareacut);
99     MCjetCont->SetRhoName(nMCrho);
100     MCjetCont->ConnectParticleContainer(MCpartCont);
101     MCjetCont->ConnectClusterContainer(MCclusCont);
102   }
103
104   AliParticleContainer *embPartCont = jetTask->AddParticleContainer(nembtracks);
105   if (embPartCont) {
106     embPartCont->SetName("EmbTracks");
107     embPartCont->SetParticlePtCut(trackptcut);
108   }
109
110   AliClusterContainer *embClusCont = jetTask->AddClusterContainer(nembclusters);
111   if (embClusCont) {
112     embClusCont->SetName("EmbClusters");
113     embClusCont->SetClusPtCut(clusptcut);
114   }
115
116   AliJetContainer *embJetCont = jetTask->AddJetContainer(nembjets,type,jetradius);
117   if (embJetCont) {
118     embJetCont->SetName("EmbJets");
119     embJetCont->SetPercAreaCut(jetareacut);
120     embJetCont->SetRhoName(nEMBrho);
121     embJetCont->ConnectParticleContainer(embPartCont);
122     embJetCont->ConnectClusterContainer(embClusCont);
123   }
124
125   /*AliParticleContainer *randPartCont = jetTask->AddParticleContainer(nrandtracks);
126   if (randPartCont) {
127     randPartCont->SetName("RandTracks");
128     randPartCont->SetParticlePtCut(trackptcut);
129   }
130
131   AliClusterContainer *randClusCont = jetTask->AddClusterContainer(nrandclusters);    
132   if (randClusCont) {
133     randClusCont->SetName("RandClusters");
134     randClusCont->SetClusPtCut(clusptcut);
135   }*/
136   
137   //-------------------------------------------------------
138   // Final settings, pass to manager and set the containers
139   //-------------------------------------------------------
140   
141   mgr->AddTask(jetTask);
142   
143   // Create containers for input/output
144   AliAnalysisDataContainer *cinput1  = mgr->GetCommonInputContainer();
145   TString contname(name);
146   contname += "_histos";
147   AliAnalysisDataContainer *coutput1 = mgr->CreateContainer(contname.Data(), 
148                                                             TList::Class(),AliAnalysisManager::kOutputContainer,
149                                                             Form("%s", AliAnalysisManager::GetCommonFileName()));
150   mgr->ConnectInput(jetTask, 0, cinput1);
151   mgr->ConnectOutput(jetTask, 1, coutput1);
152   
153   return jetTask;
154 }