]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/MakedNdeta.C
Transition PWG2/FORWARD -> PWGLF
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / MakedNdeta.C
1 /**
2  * @file   MakedNdeta.C
3  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
4  * @date   Wed Mar 23 09:41:56 2011
5  * 
6  * @brief  Run second pass analysis - make @f$ dN/d\eta@f$
7  * 
8  * @ingroup pwg2_forward_scripts_makers
9  */
10 /** 
11  * Run second pass analysis - make @f$ dN/d\eta@f$
12  * 
13  * If the ROOT AliEn interface library (libRAliEn) can be loaded, 
14  * and the parameter @a name is not empty, then use the plugin to do
15  * the analysis.  Note that in this case, the output is placed 
16  * in a sub-directory named by @a name after escaping spaces and special 
17  * characters 
18  * 
19  * @param aoddir     AOD input directory. Any file matching the pattern 
20  *                   *AliAODs*.root are added to the chain 
21  * @param nEvents    Number of events to process.  If 0 or less, then 
22  *                   all events are analysed
23  * @param trig       Trigger to use 
24  * @param useCent    Whether to use centrality or not 
25  * @param scheme     Normalisation scheme 
26  * @param vzMin      Least @f$ v_z@f$ (centimeter)
27  * @param vzMax      Largest @f$ v_z@f$ (centimeter)
28  * @param proof      If larger then 1, run in PROOF-Lite mode with this 
29  *                   many number of workers. 
30  * @param name       Name of train - free form.  This will be the name
31  *                   of the output directory if the plug-in is used 
32  * @param mcfile     Final MC corrections from this, if present
33  *
34  * @ingroup pwg2_forward_dndeta
35  */
36 void MakedNdeta(const char* aoddir   = ".", 
37                 Int_t       nEvents  = -1, 
38                 const char* trig     = "INEL",
39                 Bool_t      useCent  = false,
40                 const char* scheme   = 0,
41                 Double_t    vzMin    = -10,
42                 Double_t    vzMax    = +10,
43                 Int_t       proof    = 0,
44                 const char* name     = 0,
45                 const char* mcfile   = 0)
46 {
47   if ((name && name[0] != '\0') && gSystem->Load("libRAliEn") >= 0) {
48     gROOT->SetMacroPath(Form("%s:$(ALICE_ROOT)/PWG2/FORWARD/analysis2:"
49                              "$ALICE_ROOT/ANALYSIS/macros",
50                              gROOT->GetMacroPath()));
51     gSystem->AddIncludePath("-I${ALICE_ROOT}/include");
52     gSystem->Load("libANALYSIS");
53     gSystem->Load("libANALYSISalice");
54     gROOT->LoadMacro("TrainSetup.C+");
55     MakedNdetaTrain t(name, trig, vzMin, vzMax, scheme, useCent, false);
56     t.SetDataDir(aoddir);
57     t.SetDataSet("");
58     t.SetAllowOverwrite(true);
59     t.SetProofServer(Form("workers=%d",proof));
60     t.Run(proof > 0 ? "PROOF" : "LOCAL", "FULL", nEvents, proof > 0);
61     return;
62   }
63   // --- Libraries to load -------------------------------------------
64   gROOT->Macro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadLibs.C");
65
66   // --- Check for proof mode, and possibly upload pars --------------
67   if (proof> 0) { 
68     gROOT->LoadMacro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadPars.C");
69     LoadPars(proof);
70   }
71   
72   // --- Our data chain ----------------------------------------------
73   gROOT->LoadMacro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/MakeChain.C");
74   TChain* chain = MakeChain("AOD", aoddir,true);
75   // If 0 or less events is select, choose all 
76   if (nEvents <= 0) nEvents = chain->GetEntries();
77
78   // --- Set the macro path ------------------------------------------
79   gROOT->SetMacroPath(Form("%s:$(ALICE_ROOT)/PWG2/FORWARD/analysis2:"
80                            "$ALICE_ROOT/ANALYSIS/macros",
81                            gROOT->GetMacroPath()));
82
83   // --- Creating the manager and handlers ---------------------------
84   AliAnalysisManager *mgr  = new AliAnalysisManager(name, "Forward dN/deta");
85   AliAnalysisManager::SetCommonFileName("forward_dndeta.root");
86
87   // --- ESD input handler -------------------------------------------
88   AliAODInputHandler *aodInputHandler = new AliAODInputHandler();
89   mgr->SetInputEventHandler(aodInputHandler);      
90        
91   // --- Add tasks ---------------------------------------------------
92   // Forward 
93   gROOT->LoadMacro("AddTaskForwarddNdeta.C");
94   AddTaskForwarddNdeta(trig, vzMin, vzMax, useCent, scheme, true, mcfile);
95   // Central
96   gROOT->LoadMacro("AddTaskCentraldNdeta.C");
97   AddTaskCentraldNdeta(trig, vzMin, vzMax, useCent, scheme,false, mcfile);
98   // MC
99   gROOT->LoadMacro("AddTaskMCTruthdNdeta.C");
100   AddTaskMCTruthdNdeta(trig, vzMin, vzMax, useCent, scheme);
101
102   
103   // --- Run the analysis --------------------------------------------
104   TStopwatch t;
105   if (!mgr->InitAnalysis()) {
106     Error("MakedNdeta", "Failed to initialize analysis train!");
107     return;
108   }
109   // Skip terminate if we're so requested and not in Proof or full mode
110   mgr->SetSkipTerminate(false);
111   // Some informative output 
112   mgr->PrintStatus();
113   // mgr->SetDebugLevel(3);
114   if (mgr->GetDebugLevel() < 1 && !proof) 
115     mgr->SetUseProgressBar(kTRUE,100);
116   
117   // Run the train 
118   t.Start();
119   Printf("=== RUNNING ANALYSIS ==================================");
120   mgr->StartAnalysis(proof ? "proof" : "local", chain, nEvents);
121   t.Stop();
122   t.Print();
123 }
124 //
125 // EOF
126 //