]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/CaloTrackCorrelations/macros/AddTaskPi0.C
remove right UE histograms, duplication of default UE histograms; do the 2pi shift...
[u/mrichter/AliRoot.git] / PWGGA / CaloTrackCorrelations / macros / AddTaskPi0.C
1
2 Bool_t  kPrint         = kFALSE;
3 Bool_t  kSimulation    = kFALSE;
4 Bool_t  kUseKinematics = kFALSE;
5 Bool_t  kOutputAOD     = kFALSE;
6 Bool_t  kEventSelection= kFALSE;
7 Bool_t  kExotic        = kTRUE;
8 Bool_t  kNonLinearity  = kFALSE;
9 Int_t   kYears         = 2011;
10 TString kCollisions    = "pp";
11 TString kTrig          = "EMC7" ;
12 TString kClusterArray  = "";
13 TString kData          = ""; // MC or deltaAOD
14 TString kInputDataType = "ESD";
15 TString kCalorimeter   = "EMCAL";
16 Bool_t  kTM            = kTRUE;
17 Bool_t  kRecalTM       = kTRUE;
18 Int_t   kMinCen        = -1;
19 Int_t   kMaxCen        = -1;
20 TString kName          = "";
21 Int_t   kDebug         = -1; 
22 Float_t kVzCut         = 10.;
23 Bool_t  kPrimaryVertex = kTRUE;
24 Bool_t  kUseOADB       = kTRUE;
25 Bool_t  kTender        = kFALSE;
26
27
28 AliAnalysisTaskCaloTrackCorrelation *AddTaskPi0(const TString data          = "",
29                                                 const TString calorimeter   = "EMCAL", 
30                                                 const Bool_t  simulation    = kFALSE,
31                                                 const Bool_t  eventsel      = kFALSE,
32                                                 const Float_t vzcut         = 10,
33                                                 const Bool_t  primver       = kTRUE,
34                                                 const Bool_t  oadb          = kTRUE,
35                                                 const Bool_t  exotic        = kTRUE,
36                                                 const Bool_t  nonlin        = kFALSE,
37                                                 TString       outputfile    = "",
38                                                 const Int_t   year          = 2012,
39                                                 const TString col           = "pp", 
40                                                 const TString trigger       = "MB", 
41                                                 const TString clustersArray = "V1",
42                                                 const Bool_t  recaltm       = kTRUE,
43                                                 const Bool_t  tm            = kTRUE,
44                                                 const Int_t   minCen        = -1,
45                                                 const Int_t   maxCen        = -1,
46                                                 const Bool_t  qaan          = kFALSE,
47                                                 const Bool_t  splitan       = kFALSE,
48                                                 const Bool_t  tender        = kFALSE,
49                                                 const Bool_t  outputAOD     = kFALSE, 
50                                                 const Bool_t  printSettings = kFALSE,
51                                                 const Double_t scaleFactor   = -1
52                                                 )
53 {
54   // Creates a CaloTrackCorr task, configures it and adds it to the analysis manager.
55   
56   kPrint         = printSettings;
57   kSimulation    = simulation;
58   kYears         = year;
59   kCollisions    = col;
60   kExotic        = exotic;
61   kNonLinearity  = nonlin;
62   kTrig          = trigger;
63   kClusterArray  = clustersArray;
64   kData          = data;
65   kCalorimeter   = calorimeter;
66   kOutputAOD     = outputAOD;
67   kTM            = tm;
68   kRecalTM       = recaltm;
69   kMinCen        = minCen;
70   kMaxCen        = maxCen;
71   kEventSelection= eventsel;
72   kVzCut         = vzcut;
73   kPrimaryVertex = primver;
74   kUseOADB       = oadb;
75   kTender        = tender;
76
77   // Get the pointer to the existing analysis manager via the static access method.
78   
79   AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
80   if (!mgr) 
81   {
82     ::Error("AddTask", "No analysis manager to connect to.");
83     return NULL;
84   }  
85   
86   // Check the analysis type using the event handlers connected to the analysis manager.
87   
88   if (!mgr->GetInputEventHandler()) 
89   {
90     ::Error("AddTask", "This task requires an input event handler");
91     return NULL;
92   }
93   
94   kInputDataType = "AOD";
95   if(!kData.Contains("delta"))
96     kInputDataType = mgr->GetInputEventHandler()->GetDataType(); // can be "ESD" or "AOD"
97   
98   if(kSimulation) 
99   { 
100     kUseKinematics = (mgr->GetMCtruthEventHandler())?kTRUE:kFALSE; 
101     if (!kUseKinematics && data=="AOD" && kInputDataType != "ESD") kUseKinematics = kTRUE; //AOD primary should be available ... 
102   } 
103   
104   cout<<"********* ACCESS KINE? "<<kUseKinematics<<endl;
105   
106   // Name for containers
107   
108   kName = Form("%s_Trig%s_Cl%s_TM%d",kCalorimeter.Data(), kTrig.Data(),kClusterArray.Data(),kTM);
109
110   if(kCollisions=="PbPb" && kMaxCen>=0) kName+=Form("Cen%d_%d",kMinCen,kMaxCen);
111     
112   printf("<<<< NAME: %s >>>>>\n",kName.Data());
113   
114   // #### Configure analysis ####
115     
116   AliAnaCaloTrackCorrMaker * maker = new AliAnaCaloTrackCorrMaker();
117   
118   maker->SetScaleFactor(scaleFactor); // for MC, negative (not scaled) by default
119   
120   // General frame setting and configuration
121   maker->SetReader   (ConfigureReader()   ); 
122   maker->SetCaloUtils(ConfigureCaloUtils()); 
123   
124   // Analysis tasks setting and configuration
125   Int_t n = 0;//Analysis number, order is important  
126   
127   maker->AddAnalysis(ConfigurePhotonAnalysis(), n++); // Photon cluster selection
128   maker->AddAnalysis(ConfigurePi0Analysis(), n++); // Pi0 invariant mass accumulate 
129   maker->AddAnalysis(ConfigurePi0EbEAnalysis("Pi0", AliAnaPi0EbE::kIMCalo), n++); // Pi0 event by event selection, and photon tagging from decay    
130
131   if(qaan)                  maker->AddAnalysis(ConfigureQAAnalysis(),n++);
132   if(splitan && 
133      kCalorimeter=="EMCAL") maker->AddAnalysis(ConfigureInClusterIMAnalysis(0.5,3), n++); 
134   
135   maker->SetAnaDebug(kDebug)  ;
136   maker->SwitchOnHistogramsMaker()  ;
137   if(kData.Contains("delta")) maker->SwitchOffAODsMaker() ;
138   else                        maker->SwitchOnAODsMaker()  ;
139   
140   if(kPrint) maker->Print("");
141   
142   printf("<< End Configuration of %d analysis for calorimeter %s >>\n",n, kCalorimeter.Data());
143
144   // Create task
145   
146   AliAnalysisTaskCaloTrackCorrelation * task = new AliAnalysisTaskCaloTrackCorrelation (Form("CaloTrackCorr%s",kName.Data()));
147   task->SetConfigFileName(""); //Don't configure the analysis via configuration file.
148   task->SetDebugLevel(kDebug);
149   task->SetBranches("ESD:AliESDRun.,AliESDHeader"); 
150   task->SetAnalysisMaker(maker);
151   mgr->AddTask(task);
152   
153   //Create containers
154   
155   if(outputfile.Length()==0) outputfile = AliAnalysisManager::GetCommonFileName(); 
156   
157   AliAnalysisDataContainer *cout_pc   = mgr->CreateContainer(kName, TList::Class(), 
158                                                              AliAnalysisManager::kOutputContainer, 
159                                                              Form("%s",outputfile.Data()));
160         
161   AliAnalysisDataContainer *cout_cuts = mgr->CreateContainer(Form("Param_%s",kName.Data()), TList::Class(), 
162                                                              AliAnalysisManager::kParamContainer, 
163                                                              "AnalysisParameters.root");
164   
165   // Create ONLY the output containers for the data produced by the task.
166   // Get and connect other common input/output containers via the manager as below
167   //==============================================================================
168   mgr->ConnectInput  (task, 0, mgr->GetCommonInputContainer());
169   // AOD output slot will be used in a different way in future
170   if(!kData.Contains("delta")   && outputAOD) mgr->ConnectOutput (task, 0, mgr->GetCommonOutputContainer());
171   mgr->ConnectOutput (task, 1, cout_pc);
172   mgr->ConnectOutput (task, 2, cout_cuts);
173   
174   
175   if(kTrig=="EMC7")
176   {
177     printf("Pi0 analysis, trigger EMC7\n");
178     task->SelectCollisionCandidates(AliVEvent::kEMC7);
179   }
180   else if (kTrig=="INT7")
181   {
182     printf("Pi0 analysis, trigger INT7\n");
183     task->SelectCollisionCandidates(AliVEvent::kINT7);
184   }
185   else if(kTrig=="EMC1")
186   {
187     printf("Pi0 analysis, trigger EMC1\n");
188     task->SelectCollisionCandidates(AliVEvent::kEMC1);
189   }
190   else if(kTrig=="MB")
191   {
192     printf("Pi0 analysis, trigger MB\n");
193     task->SelectCollisionCandidates(AliVEvent::kMB);
194   }  
195   else if(kTrig=="PHOS")
196   {
197     printf("Pi0 analysis, trigger PHOS\n");
198     task->SelectCollisionCandidates(AliVEvent::kPHI7);
199   }  
200   else if(kTrig=="PHOSPb")
201   {
202     printf("Pi0 analysis, trigger PHOSPb\n");
203     task->SelectCollisionCandidates(AliVEvent::kPHOSPb);
204   }
205   else if(kTrig=="AnyINT")
206   {
207     printf("Pi0 analysis, trigger AnyINT\n");
208     task->SelectCollisionCandidates(AliVEvent::kAnyINT);
209   }  
210   else if(kTrig=="INT")
211   {
212     printf("Pi0 analysis, trigger AnyINT\n");
213     task->SelectCollisionCandidates(AliVEvent::kAny);
214   }
215   else if(kTrig=="EMCEGA")
216   {
217     printf("Pi0 analysis, trigger EMC Gamma\n");
218     task->SelectCollisionCandidates(AliVEvent::kEMCEGA);
219   } 
220   else if(kTrig=="EMCEJE")
221   {
222     printf("Pi0 analysis, trigger EMC Jet\n");
223     task->SelectCollisionCandidates(AliVEvent::kEMCEJE);
224   }
225   else if(kTrig=="Central")
226   {
227     printf("Pi0 analysis, trigger Central\n");
228     task->SelectCollisionCandidates(AliVEvent::kCentral);
229   } 
230   else if(kTrig=="SemiCentral")
231   {
232     printf("Pi0 analysis, trigger SemiCentral\n");
233     task->SelectCollisionCandidates(AliVEvent::kSemiCentral);
234   }
235   else if(kTrig=="SemiOrCentral")
236   {
237     printf("Pi0 analysis, trigger SemiCentral Or Central\n");
238     task->SelectCollisionCandidates(AliVEvent::kSemiCentral | AliVEvent::kCentral);
239   }
240   
241   
242   return task;
243 }
244
245 //____________________________________
246 AliCaloTrackReader * ConfigureReader()
247 {
248   
249   AliCaloTrackReader * reader = 0;
250   if     (kInputDataType == "ESD"&& kData=="MC" ) 
251     reader = new AliCaloTrackMCReader();
252   else if(kInputDataType=="AOD" || kData.Contains("AOD"))   
253     reader = new AliCaloTrackAODReader();
254   else if(kInputDataType=="ESD")            
255     reader = new AliCaloTrackESDReader();
256   else 
257     printf("AliCaloTrackReader::ConfigureReader() - Data combination not known kData=%s, kInputData=%s\n",kData.Data(),kInputDataType.Data());
258   
259   
260   reader->SetDebug(kDebug);//10 for lots of messages
261   
262   //Delta AOD?
263   //reader->SetDeltaAODFileName("");
264   if(kOutputAOD) reader->SwitchOnWriteDeltaAOD()  ;
265   
266   // MC settings
267   if(kUseKinematics){
268     if(kInputDataType == "ESD"){
269       reader->SwitchOnStack();          
270       reader->SwitchOffAODMCParticles(); 
271     }
272     else if(kInputDataType == "AOD"){
273       reader->SwitchOffStack();          
274       reader->SwitchOnAODMCParticles(); 
275     }
276   }  
277   
278   //------------------------
279   // Detector input filling
280   //------------------------
281   
282   //Min cluster/track E
283   reader->SetEMCALEMin(0.3); 
284   reader->SetEMCALEMax(1000); 
285   reader->SetPHOSEMin(0.3);
286   reader->SetPHOSEMax(1000);
287   reader->SetCTSPtMin(0.2);
288   reader->SetCTSPtMax(1000);
289   
290   if(!kSimulation && kCalibT) reader->SetEMCALTimeCut(-30,30); 
291   else                        reader->SetEMCALTimeCut(-1000,1000); // Open time cut
292   
293   reader->SwitchOnFiducialCut();
294   reader->GetFiducialCut()->SetSimpleCTSFiducialCut(0.8, 0, 360) ;
295
296   // Tracks
297   reader->SwitchOffCTS();
298   
299   // Calorimeter
300   
301   reader->SetEMCALClusterListName(kClusterArray);
302   if(kClusterArray == "" && !kTender) 
303   {
304     printf("**************** Standard EMCAL clusters branch analysis **************** \n");
305     reader->SwitchOnClusterRecalculation();
306     // Check in ConfigureCaloUtils that the recalibration and bad map are ON 
307   }
308   else 
309   {
310     printf("**************** Input for analysis is Clusterizer %s **************** \n", kClusterArray.Data());
311     reader->SwitchOffClusterRecalculation();
312   }  
313   
314   //if(kCalorimeter == "EMCAL") {
315     reader->SwitchOnEMCALCells();  
316     reader->SwitchOnEMCAL();
317   //}
318   //if(kCalorimeter == "PHOS") { 
319     reader->SwitchOnPHOSCells();  
320     reader->SwitchOnPHOS();
321   //}
322   
323   // for case data="deltaAOD", no need to fill the EMCAL/PHOS cluster lists
324   if(kData.Contains("delta")){
325     reader->SwitchOffEMCAL();
326     reader->SwitchOffPHOS();
327     reader->SwitchOffEMCALCells(); 
328     reader->SwitchOffPHOSCells(); 
329   }
330   
331   //-----------------
332   // Event selection
333   //-----------------
334   
335   if(!kUseKinematics && kTrig.BeginsWith("C")) reader->SetFiredTriggerClassName(kTrig); 
336   
337   // vertex event selection
338   reader->SetZvertexCut(kVzCut);               
339   if(kPrimaryVertex)reader->SwitchOnPrimaryVertexSelection(); 
340   else              reader->SwitchOffPrimaryVertexSelection(); 
341   
342   if(kEventSelection)
343   {
344     reader->SwitchOnPileUpEventRejection();   // remove pileup by default
345     reader->SwitchOnV0ANDSelection() ;        // and besides v0 AND
346   }
347   else 
348   {
349     reader->SwitchOffPileUpEventRejection();  // remove pileup by default
350     reader->SwitchOffV0ANDSelection() ;       // and besides v0 AND
351   }
352     
353   if(kCollisions=="PbPb") 
354   {
355     // Centrality
356     reader->SetCentralityClass("V0M");
357     reader->SetCentralityOpt(10);  // 10 (c= 0-10, 10-20 ...), 20  (c= 0-5, 5-10 ...) or 100 (c= 1, 2, 3 ..)
358     reader->SetCentralityBin(kMinCen,kMaxCen); // Accept all events, if not select range
359     
360     // Event plane (only used in AliAnaPi0 for the moment)
361     reader->SetEventPlaneMethod("Q");
362   }
363     
364   if(kPrint) reader->Print("");
365   
366   return reader;
367   
368 }
369
370 //_______________________________________
371 AliCalorimeterUtils* ConfigureCaloUtils()
372 {
373   
374   AliCalorimeterUtils *cu = new AliCalorimeterUtils;
375   cu->SetDebug(kDebug);
376   
377   // Remove clusters close to borders, at least max energy cell is 1 cell away 
378   cu->SetNumberOfCellsFromEMCALBorder(1);
379   cu->SetNumberOfCellsFromPHOSBorder(2);
380   
381   // Search of local maxima in cluster
382   if(kCollisions=="pp")
383   {
384     cu->SetLocalMaximaCutE(0.1);
385     cu->SetLocalMaximaCutEDiff(0.03);
386   }
387   else 
388   {
389     cu->SetLocalMaximaCutE(0.2);
390     cu->SetLocalMaximaCutEDiff(0.03);
391   }
392   
393   cu->SwitchOffClusterPlot();
394
395   if(kRecalTM) cu->SwitchOnRecalculateClusterTrackMatching(); // Done in clusterization
396   else         cu->SwitchOffRecalculateClusterTrackMatching();
397   
398   
399   //EMCAL settings
400
401   if(!kSimulation)
402     cu->SwitchOnLoadOwnEMCALGeometryMatrices();
403   
404   AliEMCALRecoUtils * recou = cu->GetEMCALRecoUtils();
405   
406   Bool_t bCalibE = kTRUE;
407   Bool_t bCalibT = kTRUE;
408   Bool_t bBadMap = kTRUE;
409   
410   if(kUseOADB)
411   {
412     cu->SwitchOnRecalibration(); 
413     cu->SwitchOnBadChannelsRemoval() ;
414     cu->SwitchOnEMCALOADB();
415   }
416   else
417   {
418     cu->SwitchOffRecalibration(); 
419     cu->SwitchOffBadChannelsRemoval() ;
420     cu->SwitchOffEMCALOADB();
421     bCalibE = kFALSE;
422     bBadMap = kFALSE;
423     bCalibT = kFALSE;
424   }
425   
426   gROOT->LoadMacro("$ALICE_ROOT/PWGGA/EMCALTasks/macros/ConfigureEMCALRecoUtils.C");
427   ConfigureEMCALRecoUtils(recou,
428                           kSimulation,                             
429                           kExotic,
430                           kNonLinearity,
431                           bCalibE, 
432                           bBadMap,
433                           bCalibT);   
434   
435
436   printf("ConfigureCaloUtils() - EMCAL Recalibration ON? %d %d\n",recou->IsRecalibrationOn(), cu->IsRecalibrationOn());
437   printf("ConfigureCaloUtils() - EMCAL BadMap        ON? %d %d\n",recou->IsBadChannelsRemovalSwitchedOn(), cu->IsBadChannelsRemovalSwitchedOn());
438   
439  
440   
441   if( kNonLinearity ) 
442   { 
443     printf("ConfigureCaloUtils() - Apply non linearity to EMCAL\n");
444     cu->SwitchOnCorrectClusterLinearity();
445   }
446   
447   if(kCalorimeter=="PHOS")
448   {
449     if (kYears <  2014) cu->SetNumberOfSuperModulesUsed(3);
450     else                cu->SetNumberOfSuperModulesUsed(4);
451   }
452   else
453   {
454     if      (kYears == 2010) cu->SetNumberOfSuperModulesUsed(4); //EMCAL first year
455     else if (kYears <  2014) cu->SetNumberOfSuperModulesUsed(10);
456     else                     cu->SetNumberOfSuperModulesUsed(20);
457   }
458   
459   // PHOS 
460   cu->SwitchOffLoadOwnPHOSGeometryMatrices();
461     
462   if(kPrint) cu->Print("");
463   
464   return cu;
465   
466 }
467
468 //_____________________________________
469 AliAnaPhoton* ConfigurePhotonAnalysis()
470 {
471   
472   AliAnaPhoton *ana = new AliAnaPhoton();
473   ana->SetDebug(kDebug); //10 for lots of messages
474   
475   // cluster selection cuts
476   
477   ana->SwitchOffFiducialCut();
478
479   ana->SetCalorimeter(kCalorimeter);
480   
481   if(kCalorimeter == "PHOS")
482   {
483     ana->SetNCellCut(2);// At least 3 cells
484     ana->SetMinPt(0.3);
485     ana->SetMinDistanceToBadChannel(2, 4, 5);
486     ana->SetTimeCut(-2000,2000); // open cut
487   }
488   else 
489   {//EMCAL
490     ana->SetNCellCut(1);// At least 2 cells
491     ana->SetMinEnergy(0.3); // avoid mip peak at E = 260 MeV
492     ana->SetMaxEnergy(1000); 
493     ana->SetTimeCut(-1000,1000); // open cut, usual time window of [425-825] ns if time recalibration is off 
494     // restrict to less than 100 ns when time calibration is on 
495     ana->SetMinDistanceToBadChannel(2, 4, 6); 
496   }
497   
498   if(kTM)
499   {
500     ana->SwitchOnTrackMatchRejection() ;
501     ana->SwitchOffTMHistoFill() ;
502   }
503   else
504   {
505     ana->SwitchOffTrackMatchRejection() ;
506     ana->SwitchOnTMHistoFill() ;
507   }
508   
509   //PID cuts (shower shape)
510   ana->SwitchOnCaloPID(); // do PID selection, unless specified in GetCaloPID, selection not based on bayesian
511   AliCaloPID* caloPID = ana->GetCaloPID();
512   //Not used in bayesian
513   
514   //EMCAL
515   caloPID->SetEMCALLambda0CutMax(0.50); // Mild cut
516   caloPID->SetEMCALLambda0CutMin(0.10);
517   
518   caloPID->SetEMCALDEtaCut(0.025);
519   caloPID->SetEMCALDPhiCut(0.030);
520   
521   //PHOS
522   caloPID->SetPHOSDispersionCut(2.5);
523   caloPID->SetPHOSRCut(2.);
524   if(kInputData=="AOD") caloPID->SetPHOSRCut(2000.); // Open cut since dX, dZ not stored
525       
526   ana->SwitchOffFillShowerShapeHistograms();  // Filled before photon shower shape selection
527   
528   // Input / output delta AOD settings
529   
530   if(!kData.Contains("delta")) 
531   {
532     ana->SetOutputAODName(Form("Photon%s",kName.Data()));
533     ana->SetOutputAODClassName("AliAODPWG4ParticleCorrelation");
534     //ana->SetOutputAODClassName("AliAODPWG4Particle"); // use if no correlation done
535   }
536   else ana->SetInputAODName(Form("Photon%s",kName.Data()));
537   
538   //Set Histograms name tag, bins and ranges
539   
540   ana->AddToHistogramsName(Form("AnaPhoton_TM%d_",kTM));
541   SetHistoRangeAndNBins(ana->GetHistogramRanges()); // see method below
542   
543   // Number of particle type MC histograms
544   ana->FillNOriginHistograms(8);
545   ana->FillNPrimaryHistograms(4);
546   
547   ConfigureMC(ana);
548   
549   if(kPrint) ana->Print("");
550   
551   return ana;
552   
553 }
554
555 //__________________________________________________________________________________________
556 AliAnaInsideClusterInvariantMass* ConfigureInClusterIMAnalysis(Float_t l0min, Float_t l0max)
557 {
558
559   AliAnaInsideClusterInvariantMass *ana = new AliAnaInsideClusterInvariantMass();
560   ana->SetDebug(kDebug); //10 for lots of messages
561   
562   // selection cuts
563   
564   ana->SetMinEnergy(5); 
565   ana->SetMaxEnergy(200.);   
566   ana->SetMinNCells(3);
567   ana->SetM02Cut(l0min,l0max);
568   ana->SetCalorimeter(kCalorimeter);
569   
570   //ana->AddToHistogramsName(Form("AnaInClusterIM_%1.2f_%1.2f_",l0min,l0max));
571   ana->AddToHistogramsName("AnaInClusterIM_");
572
573   SetHistoRangeAndNBins(ana->GetHistogramRanges()); // see method below
574   
575   AliCaloPID* caloPID = ana->GetCaloPID();
576   caloPID->SetEMCALDEtaCut(0.025);
577   caloPID->SetEMCALDPhiCut(0.030);
578   caloPID->SetClusterSplittingM02Cut(0,100); // Do the selection in the analysis class and not in the PID method to fill SS histograms
579
580   caloPID->SetPi0MassRange(0.10, 0.18);
581   caloPID->SetEtaMassRange(0.40, 0.60);
582   caloPID->SetPhotonMassRange(0.00, 0.08);  
583   
584   ConfigureMC(ana);
585   
586   if(kPrint) ana->Print("");
587   
588   return ana;
589   
590 }
591
592 //_______________________________
593 AliAnaPi0* ConfigurePi0Analysis()
594 {
595   
596   AliAnaPi0 *ana = new AliAnaPi0();
597   
598   ana->SetDebug(kDebug);//10 for lots of messages
599   
600   // Input delta AOD settings
601   ana->SetInputAODName(Form("Photon%s",kName.Data()));
602   
603   // Calorimeter settings
604   ana->SetCalorimeter(kCalorimeter);
605   
606   //settings for pp collision mixing
607   ana->SwitchOnOwnMix(); //Off when mixing done with general mixing frame
608   
609   // Cuts 
610   if(kCalorimeter=="EMCAL") ana->SetPairTimeCut(40);
611   
612   ana->SetNAsymCuts(1); // no asymmetry cut, previous studies showed small effect. 
613   //In EMCAL assymetry cut prevents combination of assymetric decays which is the main source of pi0 at high E.
614   
615   if     (kCollisions=="pp"  ) 
616   {
617     ana->SetNCentrBin(1);
618     ana->SetNZvertBin(10);
619     ana->SetNRPBin(1);
620     ana->SetNMaxEvMix(100);    
621     ana->SwitchOnSMCombinations();
622   }
623   else if(kCollisions=="PbPb") 
624   {
625     ana->SetNCentrBin(10);
626     ana->SetNZvertBin(10);
627     ana->SetNRPBin(4);
628     ana->SetNMaxEvMix(10);
629     ana->SwitchOffSMCombinations();
630   }
631
632   ana->SwitchOffMultipleCutAnalysis();
633
634   //Set Histograms name tag, bins and ranges
635   
636   ana->AddToHistogramsName(Form("AnaPi0_TM%d_",kTM));
637   SetHistoRangeAndNBins(ana->GetHistogramRanges()); // see method below
638
639   ConfigureMC(ana);
640   
641   if(kPrint) ana->Print("");
642   
643   return ana;
644   
645 }
646
647 //_____________________________________________________
648 AliAnaPi0EbE* ConfigurePi0EbEAnalysis(TString particle, 
649                                       Int_t analysis)
650 {
651   
652   AliAnaPi0EbE *ana = new AliAnaPi0EbE();
653   ana->SetDebug(kDebug);//10 for lots of messages
654   
655   ana->SetAnalysisType(analysis);
656   TString opt = "";
657   if(analysis==AliAnaPi0EbE::kIMCaloTracks) opt = "Conv";
658   if(analysis==AliAnaPi0EbE::kSSCalo)       opt = "SS";
659   
660   ana->SetMinPt(0.5);
661   
662   if(kCalorimeter=="EMCAL")ana->SetPairTimeCut(15); // More strict than in pi0 inv mass analysis
663   
664   ana->SetCalorimeter(kCalorimeter);
665   
666   // Input / output delta AOD settings
667   
668   ana->SetInputAODName(Form("Photon%s",kName.Data()));
669   if(!kInputDataType.Contains("delta")) 
670   {
671     ana->SetOutputAODName(Form("%s%s%s",particle.Data(), opt.Data(), kName.Data()));
672     ana->SetOutputAODClassName("AliAODPWG4ParticleCorrelation");
673   }
674   else  
675     ana->SetInputAODName(Form("%s%s%s",particle.Data(),opt.Data(),kName.Data()));
676   
677   if(analysis == AliAnaPi0EbE::kIMCaloTracks) ana->SetInputAODGammaConvName("PhotonsCTS");
678   
679   if(analysis!=AliAnaPi0EbE::kSSCalo)
680   {
681     AliNeutralMesonSelection *nms = ana->GetNeutralMesonSelection();
682     nms->SetParticle(particle);
683     
684     // Tighten a bit mass cut with respect to default window
685     if(particle=="Pi0") nms->SetInvMassCutRange(0.120,0.150);
686     if(particle=="Eta") nms->SetInvMassCutRange(0.520,0.580);    
687     
688     nms->SwitchOnAngleSelection();
689     nms->KeepNeutralMesonSelectionHistos(kTRUE);
690     //nms->SetAngleMaxParam(2,0.2);
691     nms->SetHistoERangeAndNBins(0, 20, 80) ;
692     //nms->SetHistoIMRangeAndNBins(0, 1, 400);
693   }
694   else
695   { // cluster splitting settings
696     ana->SetTimeCut(-1000,1000); // Open time cut
697     AliCaloPID* caloPID = ana->GetCaloPID();
698     caloPID->SetPi0MassRange(0.10, 0.18);
699     caloPID->SetEtaMassRange(0.40, 0.60);
700     caloPID->SetPhotonMassRange(0.00, 0.08);
701     caloPID->SetClusterSplittingM02Cut(0.5,100); // Do the selection in the analysis class and not in the PID method to fill SS histograms
702   }
703   
704   ana->SwitchOffSelectedClusterHistoFill(); 
705   ana->SwitchOffFillWeightHistograms();
706   
707   if(!kTM) ana->SwitchOnTMHistoFill();
708   else     ana->SwitchOffTMHistoFill();
709   
710   //Set Histograms name tag, bins and ranges
711   
712   ana->AddToHistogramsName(Form("Ana%s%sEbE_TM%d_",particle.Data(),opt.Data(),kTM));
713   SetHistoRangeAndNBins(ana->GetHistogramRanges()); // see method below
714   
715   ConfigureMC(ana);
716   
717   if(kPrint) ana->Print("");
718   
719   return  ana;
720   
721 }
722
723 //________________________________________
724 AliAnaCalorimeterQA* ConfigureQAAnalysis()
725 {
726   
727   AliAnaCalorimeterQA *ana = new AliAnaCalorimeterQA();
728   ana->SetDebug(kDebug); //10 for lots of messages
729   ana->SetCalorimeter(kCalorimeter);
730   
731   ana->SetTimeCut(-1000,1000); // Open time cut
732   
733   // Study inter detector correlation (PHOS, EMCAL, Tracks, V0)
734   if(kCalorimeter=="PHOS" && kTrig=="PHOS"){
735     ana->SwitchOnCorrelation(); // make sure you switch in the reader PHOS and EMCAL cells and clusters if option is ON
736   }
737   if(kCalorimeter=="EMCAL" && kClusterArray==""){
738     ana->SwitchOnCorrelation(); // make sure you switch in the reader PHOS and EMCAL cells and clusters if option is ON
739   }
740   else {
741     ana->SwitchOffCorrelation();
742   }
743   
744   // Study exotic clusters PHOS and EMCAL
745   if(kClusterArray==""){
746     ana->SwitchOnStudyBadClusters() ; 
747   }
748   else {
749     ana->SwitchOffStudyBadClusters() ;
750   }
751   
752   ana->SwitchOffFiducialCut();
753   ana->SwitchOffFillAllTH3Histogram();
754   ana->SwitchOffFillAllPositionHistogram();
755   ana->SwitchOffFillAllPositionHistogram2();
756   if(!kExotic)ana->SwitchOnStudyBadClusters();
757   else        ana->SwitchOffStudyBadClusters();
758   ana->SwitchOffStudyClustersAsymmetry();
759   ana->SwitchOffStudyWeight();
760   ana->SwitchOnFillAllTrackMatchingHistogram();
761   
762   ana->AddToHistogramsName("QA_"); //Begining of histograms name
763   SetHistoRangeAndNBins(ana->GetHistogramRanges()); // see method below
764   
765   ConfigureMC(ana);
766   
767   if(kPrint) ana->Print("");    
768   
769   return ana;
770   
771 }
772
773 //________________________________________________________
774 void ConfigureMC(AliAnaCaloTrackCorrBaseClass* ana)
775 {
776   if(kSimulation) ana->SwitchOnDataMC() ;//Access MC stack and fill more histograms, AOD MC not implemented yet.
777   else            ana->SwitchOffDataMC() ;
778
779   //Set here generator name, default pythia
780   //ana->GetMCAnalysisUtils()->SetMCGenerator("");
781 }  
782
783 //________________________________________________________
784 void SetHistoRangeAndNBins (AliHistogramRanges* histoRanges)
785 {
786   // Set common bins for all analysis and MC histograms filling
787     
788   histoRanges->SetHistoPtRangeAndNBins(0, 100, 200) ; // Energy and pt histograms
789   
790   if(kCalorimeter=="EMCAL")
791   {
792     if(kYears==2010)
793     {
794       histoRanges->SetHistoPhiRangeAndNBins(78*TMath::DegToRad(), 122*TMath::DegToRad(), 78) ;
795       histoRanges->SetHistoXRangeAndNBins(-230,90,120); // QA
796       histoRanges->SetHistoYRangeAndNBins(370,450,40);  // QA
797     }
798     else if(kYears==2011)
799     {           
800       histoRanges->SetHistoPhiRangeAndNBins(78*TMath::DegToRad(), 182*TMath::DegToRad(), 108) ;
801       histoRanges->SetHistoXRangeAndNBins(-600,90,200); // QA
802       histoRanges->SetHistoYRangeAndNBins(100,450,100); // QA
803     }
804     else
805     {
806       histoRanges->SetHistoPhiRangeAndNBins(78*TMath::DegToRad(), 190*TMath::DegToRad(), 122) ;
807       histoRanges->SetHistoXRangeAndNBins(-100,90,200); // QA
808       histoRanges->SetHistoYRangeAndNBins(50,450,100);  // QA
809     }
810
811     histoRanges->SetHistoEtaRangeAndNBins(-0.72, 0.72, 144) ;
812   }
813   else
814   {
815     histoRanges->SetHistoPhiRangeAndNBins(260*TMath::DegToRad(), 320*TMath::DegToRad(), 60) ;
816     histoRanges->SetHistoEtaRangeAndNBins(-0.13, 0.13, 130) ;
817   }
818   
819   histoRanges->SetHistoShowerShapeRangeAndNBins(0, 5, 500);
820   
821   // Invariant mass histoRangeslysis
822   histoRanges->SetHistoMassRangeAndNBins(0., 1., 200) ;
823   histoRanges->SetHistoAsymmetryRangeAndNBins(0., 1. , 100) ;
824   
825   // check if time calibration is on
826   histoRanges->SetHistoTimeRangeAndNBins(-1000.,1000,1000);
827   histoRanges->SetHistoDiffTimeRangeAndNBins(-200, 200, 800);
828   
829   // track-cluster residuals
830   histoRanges->SetHistoTrackResidualEtaRangeAndNBins(-0.15,0.15,300);
831   histoRanges->SetHistoTrackResidualPhiRangeAndNBins(-0.15,0.15,300);
832
833   // QA, electron, charged
834   histoRanges->SetHistoPOverERangeAndNBins(0,10.,100);
835   histoRanges->SetHistodEdxRangeAndNBins(0.,200.,200);
836   
837   // QA
838   histoRanges->SetHistoFinePtRangeAndNBins(0, 10, 200) ; // bining for fhAmpId
839   histoRanges->SetHistodRRangeAndNBins(0.,TMath::Pi(),150);
840   histoRanges->SetHistoRatioRangeAndNBins(0.,2.,100);
841   histoRanges->SetHistoVertexDistRangeAndNBins(0.,500.,500);
842   histoRanges->SetHistoNClusterCellRangeAndNBins(0,500,500);
843   histoRanges->SetHistoZRangeAndNBins(-400,400,200);
844   histoRanges->SetHistoRRangeAndNBins(400,450,25);
845   histoRanges->SetHistoV0SignalRangeAndNBins(0,5000,500);
846   histoRanges->SetHistoV0MultiplicityRangeAndNBins(0,5000,500);
847   histoRanges->SetHistoTrackMultiplicityRangeAndNBins(0,5000,500);
848   
849 }
850
851