]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGJE/EMCALJetTasks/macros/AddTaskChargedJetsPA.C
Add possibility to enable trigger patch QA in the jet preparation
[u/mrichter/AliRoot.git] / PWGJE / EMCALJetTasks / macros / AddTaskChargedJetsPA.C
1 AliAnalysisTaskChargedJetsPA* AddTaskChargedJetsPA(
2   const char*         containerSuffix         = "",
3   Double_t            jetRadius               = 0.4,
4   const char*         centralityType          = "V0A",
5   Int_t               trigger                 = AliVEvent::kINT7,
6   Bool_t              isPA                    = kTRUE,
7   Bool_t              isMC                    = kFALSE,
8   Bool_t              doJetProfileAnalysis    = kFALSE,
9   Bool_t              doTrackcutAnalysis      = kFALSE,
10   Bool_t              doJetAnalysis           = kTRUE,
11   const char*         usedTracks              = "PicoTracks",
12   Int_t               numberOfCentralityBins  = 20,
13   Double_t            areaPercentage          = 0.6,
14   Double_t            ktJetRadius             = 0.4,
15   Double_t            trackBgrdConeR          = 0.6,
16   Double_t            minJetTrackPt           = 0.150,
17   Double_t            minEta                  = -0.9,
18   Double_t            maxEta                  = +0.9,
19   Double_t            minJetEta               = -0.5,
20   Double_t            maxJetEta               = +0.5,
21   Bool_t              isEMCalTrain            = kFALSE
22 )
23 {
24   cout << " ############ MACRO EXECUTION STARTED ############\n";
25   // #### Detect the demanded trigger with its readable name
26   TString triggerName(Form("Trigger_%i", trigger));
27   if (trigger == AliVEvent::kAnyINT)
28     triggerName = "kAnyINT";
29   else if (trigger == AliVEvent::kAny)
30     triggerName = "kAny";
31   else if(trigger == AliVEvent::kINT7)
32     triggerName = "kINT7";
33   else if(trigger == AliVEvent::kMB)
34     triggerName = "kMB";
35   else if(trigger == AliVEvent::kEMC7)
36     triggerName = "kEMC7";
37   else if(trigger == AliVEvent::kEMCEJE)
38     triggerName = "kEMCEJE";
39   else if(trigger == AliVEvent::kEMCEGA)
40     triggerName = "kEMCEGA";
41
42   // #### Define manager and data container names
43   AliAnalysisManager *manager = AliAnalysisManager::GetAnalysisManager();
44   if (!manager) {
45     ::Error("AddTaskChargedJetsPA", "No analysis manager to connect to.");
46     return NULL;
47   }
48
49   TString containerNameSuffix("");
50
51   if (strcmp(containerSuffix,""))
52     containerNameSuffix = Form("_%s", containerSuffix);
53
54   TString bgrdName("Background");
55   TString myContName("");
56   TString myContJPName("");
57   TString myContTCName("");
58   if(isMC)
59   {
60     bgrdName = Form("BackgroundR0%2.0f_%s_MC%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
61     myContName = Form("AnalysisR0%2.0f_%s_MC%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
62     myContJPName = Form("JetProfileR0%2.0f_%s_MC%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
63     myContTCName = Form("TrackcutsR0%2.0f_%s_MC%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
64   }
65   else
66   {
67     bgrdName = Form("BackgroundR0%2.0f_%s%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
68     myContName = Form("AnalysisR0%2.0f_%s%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
69     myContJPName = Form("JetProfileR0%2.0f_%s%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
70     myContTCName = Form("TrackcutsR0%2.0f_%s%s", jetRadius*100, triggerName.Data(), containerNameSuffix.Data());
71   }
72
73   if(doJetAnalysis)
74   {
75     // #### Add necessary jet finder tasks
76     gROOT->LoadMacro("$ALICE_ROOT/PWGJE/EMCALJetTasks/macros/AddTaskEmcalJet.C");
77     AliEmcalJetTask* jetFinderTask = AddTaskEmcalJet(usedTracks,"",1,jetRadius,1,minJetTrackPt,0.300); // anti-kt
78     AliEmcalJetTask* jetFinderTaskKT = AddTaskEmcalJet(usedTracks,"",0,ktJetRadius,1,minJetTrackPt,0.300); // kt
79     cout << " Jet finder tasks added: " <<  jetFinderTask << " + " <<  jetFinderTaskKT << endl;
80
81     // #### Define external rho task
82     AliEmcalJetTask* jetFinderRho = AddTaskEmcalJet(usedTracks,"",1,0.4,1,minJetTrackPt,0.300); // anti-kt
83     AliEmcalJetTask* jetFinderRhoKT = AddTaskEmcalJet(usedTracks,"",0,0.4,1,minJetTrackPt,0.300); // kt
84     cout << " Jet finder tasks (used for bgrd) added: " <<  jetFinderRho << " + " <<  jetFinderRhoKT << endl;
85     gROOT->LoadMacro("$ALICE_ROOT/PWGJE/EMCALJetTasks/macros/AddTaskRhoSparse.C");
86     AliAnalysisTaskRhoSparse* rhotask = AddTaskRhoSparse(jetFinderRhoKT->GetName(), NULL, usedTracks, "", bgrdName.Data(), 0.4,"TPC", 0., 5., 0, 0,2,kFALSE,bgrdName.Data(),kTRUE);
87     cout << " Background task added: " <<  rhotask << endl;
88   }
89
90   // #### Define analysis task
91   AliAnalysisTaskChargedJetsPA *task = NULL;
92   AliAnalysisDataContainer* contHistos = manager->CreateContainer(myContName.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, Form("%s:ChargedJets", AliAnalysisManager::GetCommonFileName()));
93   if(doJetProfileAnalysis)
94     AliAnalysisDataContainer* contJetProfile = manager->CreateContainer(myContJPName.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, Form("%s:ChargedJets", AliAnalysisManager::GetCommonFileName()));
95   if(doTrackcutAnalysis)
96     AliAnalysisDataContainer* contTrackcuts = manager->CreateContainer(myContTCName.Data(), TList::Class(), AliAnalysisManager::kOutputContainer, Form("%s:ChargedJets", AliAnalysisManager::GetCommonFileName()));
97
98   if(doJetAnalysis)
99   {
100     task = new AliAnalysisTaskChargedJetsPA(Form("AnalysisPA%s_%s_%s", containerNameSuffix.Data(), jetFinderTask->GetName(), triggerName.Data()), usedTracks, jetFinderTask->GetName(),jetFinderTaskKT->GetName(), doJetProfileAnalysis, doTrackcutAnalysis);
101     task->SetExternalRhoTaskName(bgrdName.Data());
102   }
103   else
104     task = new AliAnalysisTaskChargedJetsPA(Form("AnalysisPA%s_%s", containerNameSuffix.Data(), triggerName.Data()), usedTracks, "","", doJetProfileAnalysis, doTrackcutAnalysis);
105
106   cout << " Main task created: " <<  task << endl;
107
108
109   // #### Task preferences
110   task->SetIsPA(isPA);
111   task->SetAnalyzeJetProfile(doJetProfileAnalysis);
112   task->SetAnalyzeTrackcuts(doTrackcutAnalysis);
113   task->SetAcceptanceEta(minEta,maxEta);
114   task->SetAcceptanceJetEta(minJetEta,maxJetEta);
115   task->SetSignalJetRadius(jetRadius);
116   task->SetBackgroundJetRadius(jetRadius);
117   task->SetSignalJetMinArea(areaPercentage*jetRadius*jetRadius*TMath::Pi());
118   task->SetRandConeRadius(jetRadius);
119   task->SelectCollisionCandidates(trigger);
120   task->SetCentralityType(centralityType);
121   task->SetNumberOfCentralityBins(numberOfCentralityBins);
122   task->SetDoJetAnalysis(doJetAnalysis);
123   cout << " Settings set." << endl;
124
125   // #### Add analysis task
126   manager->AddTask(task);
127   cout << " Task added to manager" << endl;
128   manager->ConnectInput(task, 0, manager->GetCommonInputContainer());
129   cout << " Input connected, common input container: " << manager->GetCommonInputContainer() << endl;
130   manager->ConnectOutput(task, 1, contHistos);
131   cout << " Output connected, contHistos: " << contHistos << endl;
132
133   if(doJetProfileAnalysis)
134   {
135     manager->ConnectOutput(task, 2, contJetProfile);
136   }
137   if(doTrackcutAnalysis && !doJetProfileAnalysis)
138   {
139     manager->ConnectOutput(task, 2, contTrackcuts);
140   }
141   else if(doTrackcutAnalysis && doJetProfileAnalysis)
142   {
143     manager->ConnectOutput(task, 3, contTrackcuts);
144   }
145
146   if(isEMCalTrain)
147     RequestMemory(task,200*1024);
148
149   cout << " ############ MACRO EXECUTION SUCCESSFUL, will return " << task << " ############\n";
150   return task;
151 }