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