]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/exa/sampleEsdAnalysis.C
Be sure to load mapping when needed
[u/mrichter/AliRoot.git] / HLT / exa / sampleEsdAnalysis.C
1 // $Id$
2 /**
3  * @file sampleEsdAnalysis.C
4  * @brief Example macro to run the AliHLTSampleESDAnalysisComponent in
5  * AliReconstruction.
6  *
7  * The component subscribes to the output of the default HLT reconstruction
8  * chain and extracts the ESD from the input.
9  *
10  * <pre>
11  * Usage: aliroot -b -q -l \
12  *     sampleEsdAnalysis.C'("file", "cdb", minEvent, maxEvent)'
13  *
14  * Examples:
15  *     sampleEsdAnalysis.C'("alien:///alice/data/2009/.../....root")' 
16  *     sampleEsdAnalysis.C'("raw.root", "local://$PWD", minEvent, MaxEvent)'
17  *
18  * Defaults
19  *     cdb="raw://"  -> take OCDB from GRID
20  *     minEvent=-1   -> no lower event selection
21  *     maxEvent=-1   -> no upper event selection
22  *
23  * </pre>
24  *
25  * The input file can be a file on the grid, indicated by the tag
26  * 'alien://' indicates. By default also the OCDB is set to the GRID.
27  * If either the file or the OCDB is taken from the GRID, the macros
28  * connects to the Grid in the beginning.
29  *
30  * As for the OCDB it is always a good idea to use the OCDB from the
31  * Grid as this will contain all the necessary objects and the latest
32  * calibration. The special URI 'raw://' is most advisable as it selects
33  * the storage automatically from the run number. Other options are e.g.
34  * - "alien://folder=/alice/data/2010/OCDB"
35  * - "local://$ALICE_ROOT/OCDB"
36  *
37  * Note: You need a valid GRID token, if you want to access files directly
38  * from the Grid, use 'alien-token-init' of your alien installation.
39  *
40  * @author Matthias.Richter@ift.uib.no
41  * @ingroup alihlt_tutorial
42  */
43 void sampleEsdAnalysis(const char *filename,
44                        const char *cdbURI,
45                        int minEvent=-1,
46                        int maxEvent=-1)
47 {
48   // connect to the GRID if we use a file or OCDB from the GRID
49   TString struri=cdbURI;
50   TString strfile=filename;
51   if (struri.BeginsWith("raw://") ||
52       strfile.Contains("://") && !strfile.Contains("local://")) {
53     TGrid::Connect("alien");
54   }
55
56   // Set the CDB storage location
57   AliCDBManager * man = AliCDBManager::Instance();
58   man->SetDefaultStorage(cdbURI);
59   if (struri.BeginsWith("local://")) {
60     // set specific storage for GRP entry
61     // search in the working directory and one level above, the latter
62     // follows the standard simulation setup like e.g. in test/ppbench
63     if (!gSystem->AccessPathName("GRP/GRP/Data")) {
64       man->SetSpecificStorage("GRP/GRP/Data", "local://$PWD");
65     } else if (!gSystem->AccessPathName("../GRP/GRP/Data")) {
66       man->SetSpecificStorage("GRP/GRP/Data", "local://$PWD/..");      
67     }
68   }
69
70   //////////////////////////////////////////////////////////////////////////////////////
71   //
72   // Reconstruction settings
73   AliReconstruction rec;
74
75   if (minEvent>=0 || maxEvent>minEvent) {
76     if (minEvent<0) minEvent=0;
77     if (maxEvent<minEvent) maxEvent=minEvent;
78     rec.SetEventRange(minEvent,maxEvent);
79   }
80   rec.SetRunReconstruction("HLT ITS TPC");
81   rec.SetWriteESDfriend(kFALSE);
82   rec.SetInput(filename);
83
84   // QA options
85   rec.SetRunQA(":") ;
86   //rec.SetQARefDefaultStorage("local://$ALICE_ROOT/QAref") ;
87
88   //////////////////////////////////////////////////////////////////////////////////////
89   //
90   // setup the HLT system
91   AliHLTSystem* pHLT=AliHLTPluginBase::GetInstance();
92
93   // define a configuration for the SampleESDAnalysis component
94   // arguments:
95   //  1) id of the configuartion, later used to refer to this configuration
96   //  2) id of the component to run
97   //  3) parents, here the configuration 'GLOBAL-esd-converter' of the libAliHLTGlobal
98   //  4) optional component arguments
99   AliHLTConfiguration esdanalysis("ESD-Analysis", "SampleESDAnalysis", "GLOBAL-esd-converter", "");
100
101   // set option for the HLT module in AliReconstruction
102   // arguments
103   //  - ignore-hltout : ignore the HLTOUT payload from the HLT DDLs
104   //  - libraries to be used as plugins
105   //  - loglevel=0x79 : Important, Warning, Error, Fatal
106   //  - chains=ESD-Analysis : chains to be run
107   rec.SetOption("HLT",
108                 "ignore-hltout " 
109                 "libAliHLTUtil.so libAliHLTGlobal.so libAliHLTSample.so "
110                 "loglevel=0x79 "
111                 "chains=ESD-Analysis "
112                 );
113
114   rec.SetRunPlaneEff(kFALSE);
115
116   // switch off cleanESD
117   rec.SetCleanESD(kFALSE);
118
119   AliLog::Flush();
120   rec.Run();
121
122 }
123
124 void sampleEsdAnalysis(const char *filename,
125                        int minEvent=-1,
126                        int maxEvent=-1)
127 {
128   sampleEsdAnalysis(filename, "raw://", minEvent, maxEvent);
129 }
130
131 void sampleEsdAnalysis()
132 {
133   cout << "sampleEsdAnalysis: Run HLT component 'SampleESDAnalyis' in AliReconstruction" << endl;
134   cout << " Usage: aliroot -b -q -l \\" << endl;
135   cout << "     sampleEsdAnalysis.C'(\"file\", \"cdb\", minEvent, maxEvent)'" << endl;
136   cout << "" << endl;
137   cout << " Examples:" << endl;
138   cout << "     sampleEsdAnalysis.C'(\"alien:///alice/data/2009/.../....root\")' " << endl;
139   cout << "     sampleEsdAnalysis.C'(\"raw.root\", \"local://$PWD\", minEvent, MaxEvent)'" << endl;
140   cout << "" << endl;
141   cout << " Defaults" << endl;
142   cout << "     cdb=\"raw://\"  -> take OCDB from GRID" << endl;
143   cout << "     minEvent=-1   -> no lower event selection" << endl;
144   cout << "     maxEvent=-1   -> no upper event selection" << endl;
145 }