]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/MakeELossFits.C
Refactoring for dN/deta tasks, more diagnostics histograms in event selector
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / MakeELossFits.C
1 /**
2  * @file 
3  * 
4  * @ingroup pwg2_forward_scripts
5  */
6 /** 
7  * Run a pass on ESD data to produce the energ loss fits 
8  * 
9  *
10  * @ingroup pwg2_forward_scripts
11  */
12 void MakeELossFits(const char* esddir, 
13                    Int_t       nEvents=1000, 
14                    Int_t       proof=0,
15                    Bool_t      mc=false)
16 {
17   // --- Libraries to load -------------------------------------------
18   gROOT->Macro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadLibs.C");
19
20   // --- Check for proof mode, and possibly upload pars --------------
21   if (proof > 0) {
22     gROOT->LoadMacro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadPars.C");
23     LoadPars(proof);
24   }
25   
26   // --- Our data chain ----------------------------------------------
27   gROOT->LoadMacro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/MakeESDChain.C");
28   TChain* chain = MakeESDChain(esddir, true);
29   // If 0 or less events is select, choose all 
30   if (nEvents <= 0) nEvents = chain->GetEntries();
31   Info("MakeELossFits", "Will analyse %d events", nEvents);
32
33   // --- Creating the manager and handlers ---------------------------
34   AliAnalysisManager *mgr  = new AliAnalysisManager("Forward ELoss Train", 
35                                                     "Forward energy loss");
36   AliAnalysisManager::SetCommonFileName("forward_eloss.root");
37
38   AliESDInputHandler *esdHandler = new AliESDInputHandler();
39   esdHandler->SetInactiveBranches("AliESDACORDE "
40                                   "AliRawDataErrorLogs "
41                                   "CaloClusters "
42                                   "Cascades "
43                                   "EMCALCells "
44                                   "EMCALTrigger "
45                                   "Kinks "
46                                   "Cascades "
47                                   "MuonTracks "
48                                   "TrdTracks "
49                                   "CaloClusters "
50                                   "HLTGlobalTrigger");
51   mgr->SetInputEventHandler(esdHandler);      
52   
53   // AOD output handler
54   AliAODHandler* aodHandler   = new AliAODHandler();
55   mgr->SetOutputEventHandler(aodHandler);
56   aodHandler->SetOutputFileName("AliAODs.root");
57
58   // --- Add tasks ---------------------------------------------------
59   gROOT->LoadMacro("$ALICE_ROOT/ANALYSIS/macros/AddTaskPhysicsSelection.C");
60   AddTaskPhysicsSelection(mc, kTRUE, kTRUE);
61   AliFMDEnergyFitterTask* task = new AliFMDEnergyFitterTask("fmdEnergyFitter");
62   mgr->AddTask(task);
63
64   // --- Make the output container and connect it --------------------
65   TString outputfile = AliAnalysisManager::GetCommonFileName();
66   AliAnalysisDataContainer* histOut = 
67     mgr->CreateContainer("Forward", TList::Class(), 
68                          AliAnalysisManager::kOutputContainer,outputfile);
69   mgr->ConnectInput(task, 0, mgr->GetCommonInputContainer());
70   mgr->ConnectOutput(task, 1, histOut);
71
72   // --- Set parameters on the algorithms ----------------------------
73   // Set the number of SPD tracklets for which we consider the event a
74   // low flux event
75   task->GetEventInspector().SetLowFluxCut(1000); 
76   // Set the maximum error on v_z [cm]
77   task->GetEventInspector().SetMaxVzErr(0.2);
78   // Set the eta axis to use - note, this overrides whatever is used
79   // by the rest of the algorithms - but only for the energy fitter
80   // algorithm. 
81   task->GetEnergyFitter().SetEtaAxis(200, -4, 6);
82   // Set maximum energy loss to consider 
83   task->GetEnergyFitter().SetMaxE(15); 
84   // Set number of energy loss bins 
85   task->GetEnergyFitter().SetNEbins(100);
86   // Set whether to use increasing bin sizes 
87   task->GetEnergyFitter().SetUseIncreasingBins(true);
88   // Set whether to do fit the energy distributions 
89   task->GetEnergyFitter().SetDoFits(kTRUE);
90   // Set whether to make the correction object 
91   task->GetEnergyFitter().SetDoMakeObject(kTRUE);
92   // Set the low cut used for energy
93   task->GetEnergyFitter().SetLowCut(0.4);
94   // Set the number of bins to subtract from maximum of distributions
95   // to get the lower bound of the fit range
96   task->GetEnergyFitter().SetFitRangeBinWidth(4);
97   // Set the maximum number of landaus to try to fit (max 5)
98   task->GetEnergyFitter().SetNParticles(5);
99   // Set the minimum number of entries in the distribution before
100   // trying to fit to the data
101   task->GetEnergyFitter().SetMinEntries(1000);
102   // --- Set limits on fits the energy -------------------------------
103   // Maximum relative error on parameters 
104   AliFMDCorrELossFit::ELossFit::fgMaxRelError = .12;
105   // Least weight to use 
106   AliFMDCorrELossFit::ELossFit::fgLeastWeight = 1e-5;
107   // Maximum value of reduced chi^2 
108   AliFMDCorrELossFit::ELossFit::fgMaxChi2nu   = 10;
109   
110   // --- Run the analysis --------------------------------------------
111   TStopwatch t;
112   if (!mgr->InitAnalysis()) {
113     Error("RunManager", "Failed to initialize analysis train!");
114     return;
115   }
116   // Some informative output 
117   mgr->PrintStatus();
118   // mgr->SetDebugLevel(3);
119   if (mgr->GetDebugLevel() < 1 && proof <= 0) mgr->SetUseProgressBar(kTRUE);
120
121   // Run the train 
122   t.Start();
123   Printf("=== RUNNING ANALYSIS ==================================");
124   mgr->StartAnalysis(proof > 0 ? "proof" : "local", chain, nEvents);
125   t.Stop();
126   t.Print();
127 }
128 //
129 // EOF
130 //