]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/macros/AddTaskPi0.C
Many fixes:
[u/mrichter/AliRoot.git] / PWG4 / macros / AddTaskPi0.C
1 AliAnalysisTaskParticleCorrelation *AddTaskPi0(TString data, TString calorimeter, Bool_t kPrintSettings = kFALSE,Bool_t kSimulation = kFALSE, Bool_t outputAOD=kFALSE, Bool_t oldAOD=kFALSE)
2 {
3   // Creates a PartCorr task, configures it and adds it to the analysis manager.
4   
5   // Get the pointer to the existing analysis manager via the static access method.
6   //==============================================================================
7   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
8   if (!mgr) {
9     ::Error("AddTaskPi0", "No analysis manager to connect to.");
10     return NULL;
11   }  
12   
13   // Check the analysis type using the event handlers connected to the analysis manager.
14   //==============================================================================
15   if (!mgr->GetInputEventHandler()) {
16     ::Error("AddTaskPi0", "This task requires an input event handler");
17     return NULL;
18   }
19   TString inputDataType = "AOD";
20   if(!data.Contains("delta"))
21     inputDataType = mgr->GetInputEventHandler()->GetDataType(); // can be "ESD" or "AOD"
22   //cout<<"DATA TYPE :: "<<inputDataType<<endl;
23   // inputDataType: data managed by the input handler
24   // data: can be same as one managed by input handler, or the output AOD created by the filter. By default use AOD
25   
26   Bool_t kUseKinematics = kFALSE; 
27   if(kSimulation) { 
28     kUseKinematics = (mgr->GetMCtruthEventHandler())?kTRUE:kFALSE; 
29     if (!kUseKinematics && data=="AOD" && inputDataType != "ESD") kUseKinematics = kTRUE; //AOD primary should be available ... 
30   } 
31   
32   cout<<"********* ACCESS KINE? "<<kUseKinematics<<endl;
33   
34   // Configure analysis
35   //===========================================================================
36   
37   // *** Reader ***
38   AliCaloTrackReader * reader = ;
39   if(data.Contains("AOD")) reader = new AliCaloTrackAODReader();
40   else if(data=="ESD") reader = new AliCaloTrackESDReader();
41   else if(data=="MC" && inputDataType == "ESD") reader = new AliCaloTrackMCReader();
42   reader->SetDebug(-1);//10 for lots of messages
43   reader->SwitchOnCTS();
44   //reader->SetDeltaAODFileName("");
45   //if(!kSimulation) reader->SetFiredTriggerClassName("CINT1B-ABCE-NOPF-ALL");
46   if(calorimeter == "EMCAL") {
47     reader->SwitchOnEMCALCells();  
48     reader->SwitchOnEMCAL();
49   }
50   if(calorimeter == "PHOS") { 
51     reader->SwitchOnPHOSCells();  
52     reader->SwitchOnPHOS();
53   }
54   
55   // for case data="deltaAOD", no need to fill the EMCAL/PHOS cluster lists
56   if(data.Contains("delta")){
57     reader->SwitchOffEMCAL();
58     reader->SwitchOffPHOS();
59     reader->SwitchOffEMCALCells(); 
60     reader->SwitchOffPHOSCells(); 
61   }
62   
63   if(kUseKinematics){
64     if(inputDataType == "ESD"){
65       reader->SwitchOnStack();          
66       reader->SwitchOffAODMCParticles(); 
67     }
68     else if(inputDataType == "AOD"){
69       reader->SwitchOffStack();          
70       reader->SwitchOnAODMCParticles(); 
71     }
72   }
73   
74   //Min particle pT
75   reader->SetEMCALPtMin(0.1); 
76   reader->SetPHOSPtMin(0.);
77   reader->SetCTSPtMin(0.);
78   if(outputAOD)  reader->SwitchOnWriteDeltaAOD()  ;
79   if(oldAOD) reader->SwitchOnOldAODs();
80   if(kPrintSettings) reader->Print("");
81   
82   // *** Calorimeters Utils     ***
83   AliCalorimeterUtils *cu = new AliCalorimeterUtils;
84   // Remove clusters close to borders, at least max energy cell is 1 cell away 
85   cu->SetNumberOfCellsFromEMCALBorder(1);
86   cu->SetNumberOfCellsFromPHOSBorder(2);
87   
88   // Remove EMCAL hottest channels for first LHC10 periods      
89   cu->SwitchOnBadChannelsRemoval();
90   // SM0
91   cu->SetEMCALChannelStatus(0,3,13);  cu->SetEMCALChannelStatus(0,44,1); cu->SetEMCALChannelStatus(0,3,13); 
92   cu->SetEMCALChannelStatus(0,20,7);  cu->SetEMCALChannelStatus(0,38,2);   
93   // SM1
94   cu->SetEMCALChannelStatus(1,4,7);   cu->SetEMCALChannelStatus(1,4,13);  cu->SetEMCALChannelStatus(1,9,20); 
95   cu->SetEMCALChannelStatus(1,14,15); cu->SetEMCALChannelStatus(1,23,16); cu->SetEMCALChannelStatus(1,32,23); 
96   cu->SetEMCALChannelStatus(1,37,5);  cu->SetEMCALChannelStatus(1,40,1);  cu->SetEMCALChannelStatus(1,40,2);
97   cu->SetEMCALChannelStatus(1,40,5);  cu->SetEMCALChannelStatus(1,41,0);  cu->SetEMCALChannelStatus(1,41,1);
98   cu->SetEMCALChannelStatus(1,41,2);  cu->SetEMCALChannelStatus(1,41,4);
99   // SM2        
100   cu->SetEMCALChannelStatus(2,14,15); cu->SetEMCALChannelStatus(2,18,16); cu->SetEMCALChannelStatus(2,18,17); 
101   cu->SetEMCALChannelStatus(2,18,18); cu->SetEMCALChannelStatus(2,18,20); cu->SetEMCALChannelStatus(2,18,21); 
102   cu->SetEMCALChannelStatus(2,18,23); cu->SetEMCALChannelStatus(2,19,16); cu->SetEMCALChannelStatus(2,19,17); 
103   cu->SetEMCALChannelStatus(2,19,19); cu->SetEMCALChannelStatus(2,19,20); cu->SetEMCALChannelStatus(2,19,21); 
104   cu->SetEMCALChannelStatus(2,19,22);
105   //SM3
106   cu->SetEMCALChannelStatus(3,4,7);
107   
108   
109   //Recalibration
110   //cu->SwitchOnRecalibration();
111   //TFile * f = new TFile("RecalibrationFactors.root","read");
112   //cu->SetEMCALChannelRecalibrationFactors(0,(TH2F*)f->Get("EMCALRecalFactors_SM0"));
113   //cu->SetEMCALChannelRecalibrationFactors(1,(TH2F*)f->Get("EMCALRecalFactors_SM1"));
114   //cu->SetEMCALChannelRecalibrationFactors(2,(TH2F*)f->Get("EMCALRecalFactors_SM2"));
115   //cu->SetEMCALChannelRecalibrationFactors(3,(TH2F*)f->Get("EMCALRecalFactors_SM3"));
116   //f->Close(); 
117   
118   cu->SetDebug(-1);
119   if(kPrintSettings) cu->Print("");
120   
121   
122   // ##### Analysis algorithm settings ####
123   
124   // -------------------------------------------------
125   // --- Photon/Pi0/Omega/Electron Analysis ---
126   // -------------------------------------------------
127   
128   AliAnaPhoton *anaphoton = new AliAnaPhoton();
129   anaphoton->SetDebug(-1); //10 for lots of messages
130   if(calorimeter == "PHOS"){
131     anaphoton->SetNCellCut(0);// At least 2 cells
132     anaphoton->SetMinPt(0.);
133     anaphoton->SetMinDistanceToBadChannel(2, 4, 5);
134   }
135   else {//EMCAL
136     //anaphoton->SetNCellCut(0);// At least 2 cells
137     anaphoton->SetMinPt(0.1); // no effect minium EMCAL cut.
138     if(!kUseKinematics) anaphoton->SetTimeCut(400,900);// Time window of [400-900] ns
139     anaphoton->SetMinDistanceToBadChannel(6, 12, 18);
140   }
141   anaphoton->SetCalorimeter(calorimeter);
142   if(kUseKinematics) anaphoton->SwitchOnDataMC() ;//Access MC stack and fill more histograms
143   else  anaphoton->SwitchOffDataMC() ;
144   anaphoton->SwitchOffCaloPID();
145   anaphoton->SwitchOffFiducialCut();
146   if(kSimulation){
147     anaphoton->SwitchOnFiducialCut();
148     AliFiducialCut * fidCut1stYear = anaphoton->GetFiducialCut();
149     fidCut1stYear->DoCTSFiducialCut(kFALSE) ;
150     fidCut1stYear->DoEMCALFiducialCut(kTRUE) ;
151     fidCut1stYear->DoPHOSFiducialCut(kTRUE) ;
152     fidCut1stYear->SetSimpleEMCALFiducialCut(0.7,80.,120.);
153     fidCut1stYear->SetSimplePHOSFiducialCut(0.12,260.,320.);
154   }
155   
156   if(!data.Contains("delta")) {
157     anaphoton->SetOutputAODName(Form("Photons%s",calorimeter.Data()));
158     anaphoton->SetOutputAODClassName("AliAODPWG4ParticleCorrelation");
159   }
160   else anaphoton->SetInputAODName(Form("Photons%s",calorimeter.Data()));
161   anaphoton->AddToHistogramsName("AnaPhotonCorr_");
162   //Set Histograms bins and ranges
163   anaphoton->SetHistoPtRangeAndNBins(0, 50, 200) ;
164   //      ana->SetHistoPhiRangeAndNBins(0, TMath::TwoPi(), 100) ;
165   //      ana->SetHistoEtaRangeAndNBins(-0.7, 0.7, 100) ;
166   if(kPrintSettings) anaphoton->Print("");
167   
168   // -----------------------------------
169   // --- Pi0 Invariant Mass Analysis ---
170   // -----------------------------------
171   
172   AliAnaPi0 *anapi0 = new AliAnaPi0();
173   anapi0->SetDebug(-1);//10 for lots of messages
174   anapi0->SetInputAODName(Form("Photons%s",calorimeter.Data()));
175   anapi0->SetCalorimeter(calorimeter);
176   
177   anapi0->SwitchOnMultipleCutAnalysis(); 
178   anapi0->SetNPtCuts(2);
179   anapi0->SetNAsymCuts(2);
180   anapi0->SetNNCellCuts(2);
181   anapi0->SetNPIDBits(2);
182   
183   anapi0->SetPtCutsAt(0,0.3); anapi0->SetPtCutsAt(1,0.5);
184   anapi0->SetAsymCutsAt(0,0.1);anapi0->SetAsymCutsAt(1,0.5);
185   anapi0->SetNCellCutsAt(0,1); anapi0->SetNCellCutsAt(1,2);
186   anapi0->SetPIDBitsAt(0,2);  anapi0->SetPIDBitsAt(1,4);
187
188   
189   if(kSimulation){
190     anapi0->SwitchOnFiducialCut();
191     AliFiducialCut * fidCut1stYear = anapi0->GetFiducialCut();
192     fidCut1stYear->DoCTSFiducialCut(kFALSE) ;
193     fidCut1stYear->DoEMCALFiducialCut(kTRUE) ;
194     fidCut1stYear->DoPHOSFiducialCut(kTRUE) ;
195     fidCut1stYear->SetSimpleEMCALFiducialCut(0.7,80.,120.);
196     fidCut1stYear->SetSimplePHOSFiducialCut(0.12,260.,320.);
197   }  
198         
199   anapi0->SetNPID(1); //Available from tag AliRoot::v4-18-15-AN
200   //settings for pp collision mixing
201   anapi0->SwitchOnOwnMix(); //Off when mixing done with general mixing frame
202   anapi0->SetNCentrBin(1);
203   anapi0->SetNZvertBin(1);
204   anapi0->SetNRPBin(1);
205   anapi0->SetNMaxEvMix(10);
206   
207   if(kUseKinematics)anapi0->SwitchOnDataMC() ;//Access MC stack and fill more histograms
208   else              anapi0->SwitchOffDataMC() ;
209   if(calorimeter=="PHOS") anapi0->SetNumberOfModules(3); //PHOS first year
210   else  anapi0->SetNumberOfModules(4); //EMCAL first year
211   anapi0->SetHistoPtRangeAndNBins(0, 50, 200) ;
212   //anapi0->SetHistoPhiRangeAndNBins(0, TMath::TwoPi(), 100) ;
213   //anapi0->SetHistoEtaRangeAndNBins(-0.8, 0.8, 200) ;
214   anapi0->SetHistoMassRangeAndNBins(0., 0.6, 200) ;
215   anapi0->SetHistoAsymmetryRangeAndNBins(0., 1. , 10) ;
216   if(kPrintSettings) anapi0->Print("");
217
218   
219   
220   // #### Configure Maker ####
221   AliAnaPartCorrMaker * maker = new AliAnaPartCorrMaker();
222   maker->SetReader(reader);//pointer to reader
223   maker->SetCaloUtils(cu); //pointer to calorimeter utils
224   Int_t n = 0;//Analysis number, order is important
225   // Particle selection analysis
226   maker->AddAnalysis(anaphoton,n++);
227   maker->AddAnalysis(anapi0,n++);
228   maker->SetAnaDebug(-1)  ;
229   maker->SwitchOnHistogramsMaker()  ;
230   if(data.Contains("delta")) maker->SwitchOffAODsMaker()  ;
231   else                       maker->SwitchOnAODsMaker()  ;
232         
233   if(kPrintSettings) maker->Print("");
234   
235   printf("======================== \n");
236   printf(" End Configuration of Pi0 analysis with detector %s \n",calorimeter.Data());
237   printf("======================== \n");
238   
239   // Create task
240   //===========================================================================
241   AliAnalysisTaskParticleCorrelation * task = new AliAnalysisTaskParticleCorrelation (Form("PartCorr%s",calorimeter.Data()));
242   task->SetConfigFileName(""); //Don't configure the analysis via configuration file.
243   //task->SetDebugLevel(-1);
244   task->SelectCollisionCandidates();
245   task->SetAnalysisMaker(maker);
246   //if(!kSimulation)task->SelectCollisionCandidates(); //AliPhysicsSelection has to be attached before.
247   mgr->AddTask(task);
248   
249   //Create containers
250   char name[128];
251   sprintf(name,"PartCorr_%s",calorimeter.Data());
252   cout<<"Name of task "<<name<<endl;
253   //AliAnalysisDataContainer *cout_pc = mgr->CreateContainer(Form(name),TList::Class(),
254   //                                       AliAnalysisManager::kOutputContainer, Form("PartCorr_%s.root",calorimeter.Data()));
255   
256   TString outputfile = AliAnalysisManager::GetCommonFileName(); 
257   //  AliAnalysisDataContainer *cout_pc = mgr->CreateContainer(Form("PartCorr_%s",calorimeter.Data()),  TList::Class(), AliAnalysisManager::kOutputContainer, Form("%s:PartCorr_%s",outputfile.Data(),calorimeter.Data()));
258   AliAnalysisDataContainer *cout_pc   = mgr->CreateContainer(calorimeter.Data(), TList::Class(), 
259                                                              AliAnalysisManager::kOutputContainer, 
260                                                              Form("%s:PartCorr",outputfile.Data()));
261         
262   AliAnalysisDataContainer *cout_cuts = mgr->CreateContainer(Form("%sCuts",calorimeter.Data()), TList::Class(), 
263                                                              AliAnalysisManager::kParamContainer, 
264                                                              Form("%s:PartCorrCuts",outputfile.Data()));
265         
266   // Create ONLY the output containers for the data produced by the task.
267   // Get and connect other common input/output containers via the manager as below
268   //==============================================================================
269   mgr->ConnectInput  (task, 0, mgr->GetCommonInputContainer());
270   // AOD output slot will be used in a different way in future
271   if(!data.Contains("delta")   && outputAOD) mgr->ConnectOutput (task, 0, mgr->GetCommonOutputContainer());
272   mgr->ConnectOutput (task, 1, cout_pc);
273   mgr->ConnectOutput (task, 2, cout_cuts);
274   
275   return task;
276 }
277
278