]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/MakeMCCorr.C
Updates to scripts. Mostly documentation and some new functionalities
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / MakeMCCorr.C
1 /**
2  * @file   MakeMCCorr.C
3  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
4  * @date   Tue Jul 12 10:06:07 2011
5  * 
6  * @brief  Generate MC corrections 
7  * 
8  * @ingroup pwg2_forward_scripts_makers
9  */
10 //====================================================================
11 /** 
12  * Generatew MC corrections 
13  * 
14  * If the ROOT AliEn interface library (libRAliEn) can be loaded, 
15  * and the parameter @a name is not empty, then use the plugin to do
16  * the analysis.  Note that in this case, the output is placed 
17  * in a sub-directory named by @a name after escaping spaces and special 
18  * characters 
19  * 
20  * @param esddir     AOD input directory. Any file matching the pattern 
21  *                   *AliAODs*.root are added to the chain 
22  * @param nEvents    Number of events to process.  If 0 or less, then 
23  *                   all events are analysed
24  * @param vzMin      Least @f$ v_z@f$ (centimeter)
25  * @param vzMax      Largest @f$ v_z@f$ (centimeter)
26  * @param proof      If larger then 1, run in PROOF-Lite mode with this 
27  *                   many number of workers. 
28  * @param name       Name of train - free form.  This will be the name
29  *                   of the output directory if the plug-in is used 
30  *
31  * @ingroup pwg2_forward_corr
32  * @ingroup pwg2_forward_scripts_makers
33  */
34 void MakeMCCorr(const char* esddir   = ".", 
35                 Int_t       nEvents  = -1, 
36                 Double_t    vzMin    = -10,
37                 Double_t    vzMax    = +10,
38                 Int_t       proof    = 0,
39                 const char* name     = 0)
40 {
41   if ((name && name[0] != '\0') && gSystem->Load("libRAliEn") >= 0) {
42     Error("MakeMCCorr", "Plug-in mode not implemented yet!");
43     return;
44
45     const char* fwdPath = 
46       gSystem->ExpandPathName("$(ALICE_ROOT)/PWG2/FORWARD/analysis2");
47     gROOT->LoadMacro(Form("%s/trains/BuildTrain.C", fwdPath));
48     BuildTrain("MakeMCCorrTrain");
49
50     MakeMCCorrTrain t(name, vzMin, vzMax);
51     t.SetDataDir(esddir);
52     t.SetDataSet("");
53     t.SetAllowOverwrite(true);
54     t.SetProofServer(Form("workers=%d",proof));
55     t.Run(proof > 0 ? "PROOF" : "LOCAL", "FULL", nEvents, proof > 0);
56     return;
57   }
58   // --- Libraries to load -------------------------------------------
59   gROOT->Macro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadLibs.C");
60
61   // --- Check for proof mode, and possibly upload pars --------------
62   if (proof> 0) { 
63     gROOT->LoadMacro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadPars.C");
64     LoadPars(proof);
65   }
66   
67   // --- Our data chain ----------------------------------------------
68   gROOT->LoadMacro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/MakeChain.C");
69   TChain* chain = MakeChain("ESD", esddir,true);
70   // If 0 or less events is select, choose all 
71   if (nEvents <= 0) nEvents = chain->GetEntries();
72
73   // --- Set the macro path ------------------------------------------
74   gROOT->SetMacroPath(Form("%s:$(ALICE_ROOT)/PWG2/FORWARD/analysis2:"
75                            "$ALICE_ROOT/ANALYSIS/macros",
76                            gROOT->GetMacroPath()));
77
78   // --- Creating the manager and handlers ---------------------------
79   AliAnalysisManager *mgr  = new AliAnalysisManager(name, "Forward MC corr");
80   AliAnalysisManager::SetCommonFileName("forward_mccorr.root");
81
82   // --- ESD input handler -------------------------------------------
83   AliESDInputHandler *esdHandler = new AliESDInputHandler();
84   mgr->SetInputEventHandler(esdHandler);      
85        
86   // --- Monte Carlo handler -----------------------------------------
87   AliMCEventHandler* mcHandler = new AliMCEventHandler();
88   mgr->SetMCtruthEventHandler(mcHandler);
89   mcHandler->SetReadTR(true);    
90
91   // --- Add Physics Selection ---------------------------------------
92   // Physics selection 
93   gROOT->LoadMacro("AddTaskPhysicsSelection.C");
94   AddTaskPhysicsSelection(kTRUE, kTRUE, kFALSE);
95   // --- Fix up physics selection to give proper A,C, and E triggers -
96   AliInputEventHandler* ih =
97     static_cast<AliInputEventHandler*>(mgr->GetInputEventHandler());
98   AliPhysicsSelection* ps = 
99     static_cast<AliPhysicsSelection*>(ih->GetEventSelection());
100   // Ignore trigger class when selecting events.  This mean that we
101   // get offline+(A,C,E) events too
102   // ps->SetSkipTriggerClassSelection(true);
103
104   // --- Add our task ------------------------------------------------
105   gROOT->LoadMacro("AddTaskForwardMCCorr.C");
106   AddTaskForwardMCCorr();
107
108   gROOT->LoadMacro("AddTaskCentralMCCorr.C");
109   AddTaskCentralMCCorr();
110   
111   // --- Run the analysis --------------------------------------------
112   TStopwatch t;
113   if (!mgr->InitAnalysis()) {
114     Error("MakedNdeta", "Failed to initialize analysis train!");
115     return;
116   }
117   // Skip terminate if we're so requested and not in Proof or full mode
118   mgr->SetSkipTerminate(false);
119   // Some informative output 
120   mgr->PrintStatus();
121   // mgr->SetDebugLevel(3);
122   if (mgr->GetDebugLevel() < 1 && !proof) 
123     mgr->SetUseProgressBar(kTRUE,100);
124   
125   // Run the train 
126   t.Start();
127   Printf("=== RUNNING ANALYSIS ==================================");
128   mgr->StartAnalysis(proof ? "proof" : "local", chain, nEvents);
129   t.Stop();
130   t.Print();
131 }
132 //
133 // EOF
134 //