]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/RunLinkToMCAnalysisExample.C
c4eac3eb4ab5bc855160bd69a7264f0f80951281
[u/mrichter/AliRoot.git] / PWG3 / muon / RunLinkToMCAnalysisExample.C
1 /**************************************************************************
2  * This file is property of and copyright by the ALICE HLT Project        *
3  * All rights reserved.                                                   *
4  *                                                                        *
5  * Primary Authors:                                                       *
6  *   Artur Szostak <artursz@iafrica.com>                                  *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /// \ingroup macros
18 /// \file RunLinkToMCAnalysisExample.C
19 /// \brief Example macro for running the AliAnalysisTaskLinkToMC analysis task.
20 /// \author Artur Szostak <artursz@iafrica.com>
21 ///
22 /// This macro shows an example of how to run the AliAnalysisTaskLinkToMC analysis
23 /// task in local mode. It can be used as a quick check of the analysis task.
24 /// It will generate AliAOD.root and hists.root files. The hists.root can then
25 /// be used in the PlotEfficiency.C macro.
26 /// Run this macro as follows:
27 ///
28 ///  $ aliroot RunLinkToMCAnalysisExample.C\(\"AliESDs.root\"\)
29 ///
30 /// where AliESDs.root should be the correct path to a root file containing ESD
31 /// objects created by the offline reconstruction.
32
33
34 #if !defined(__CINT__) || defined(__MAKECINT__)
35 #error This macro must be run in interpreted mode.
36 #endif
37
38
39 void RunLinkToMCAnalysisExample(const char* esdFile = "./AliESDs.root")
40 {
41         // Load needed libraries
42         gSystem->Load("libTree.so");
43         gSystem->Load("libGeom.so");
44         gSystem->Load("libVMC.so");
45         gSystem->Load("libPhysics.so");
46         gSystem->Load("libSTEERBase.so");
47         gSystem->Load("libESD.so");
48         gSystem->Load("libAOD.so");
49         gSystem->Load("libANALYSIS.so");
50         gSystem->Load("libANALYSISalice.so");
51         gSystem->Load("libPWG3base.so");
52         gSystem->Load("libPWG3muon.so");
53         
54         // Create the TChain for esdTrees in the AliESDs.root file.
55         TChain* chain = new TChain("esdTree");
56         chain->Add(esdFile);
57         if (!chain) return;
58         
59         // Create the analysis manager and event handlers.
60         AliAnalysisManager* mgr = new AliAnalysisManager("Analysis Train", "An example analysis train setup for AliAnalysisTaskLinkToMC.");
61         AliESDInputHandler* esdHandler = new AliESDInputHandler();
62         mgr->SetInputEventHandler(esdHandler);
63         AliMCEventHandler* mcHandler = new AliMCEventHandler();
64         mgr->SetMCtruthEventHandler(mcHandler);
65         mcHandler->SetReadTR(kTRUE); 
66         AliAODHandler* aodHandler = new AliAODHandler();
67         mgr->SetOutputEventHandler(aodHandler);
68         aodHandler->SetOutputFileName("AliAOD.root");
69         
70         // Create the analysis task and setup the parameters.
71         AliAnalysisTaskLinkToMC* linktask = new AliAnalysisTaskLinkToMC("Task to link ESD tracks to corresponding MC tracks.");
72         linktask->MinClusters(6);
73         linktask->HardCutLimitX(4);
74         linktask->HardCutLimitY(4);
75         linktask->SigmaCut(5.);
76         linktask->MinClustersInSt45(3);
77         linktask->StationMustMatch(1, true);  // At least one cluster in station 1 must match.
78         linktask->StationMustMatch(2, true);  // At least one cluster in station 2 must match.
79         linktask->StationMustMatch(3, true);  // At least one cluster in station 3 must match.
80         linktask->GenerateHistograms(true);
81         mgr->AddTask(linktask);
82         
83         // Create the input and output containers and connect them up to the analysis task.
84         AliAnalysisDataContainer* cinEsd = mgr->CreateContainer("cESD", TChain::Class(), AliAnalysisManager::kInputContainer);
85         AliAnalysisDataContainer* coutAod = mgr->CreateContainer("cAOD", TTree::Class(), AliAnalysisManager::kOutputContainer, "default");
86         AliAnalysisDataContainer* coutHists = mgr->CreateContainer("cHists", TList::Class(), AliAnalysisManager::kOutputContainer, "hists.root");
87         mgr->ConnectInput(linktask, 0, cinEsd);
88         mgr->ConnectOutput(linktask, 0, coutAod);
89         mgr->ConnectOutput(linktask, 1, coutHists);
90         
91         if (mgr->InitAnalysis())
92         {
93                 mgr->PrintStatus();
94                 mgr->StartAnalysis("local", chain);
95         }
96 }
97