]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/TRD/macros/rec-hlt-trd.cxx
updating and adding macros for HLT TRD evaluation (Theodor)
[u/mrichter/AliRoot.git] / HLT / TRD / macros / rec-hlt-trd.cxx
1 // This macro is used to simulate the HLT reconstruction
2 // usage: aliroot rec-hlt-trd.cxx("/data/run/raw.root")    reconstruct local raw root file (you might add "alien://" to reconstruct remotely)
3 //    or: aliroot rec-hlt-trd.cxx("/data/run/")            reconstruct local raw ddls *1
4 //    or copy into folder and aliroot rec-hlt-trd.cxx      reconstruct raw.root in pwd
5 //
6 // (*1) here /data/run/ must contain subfolders rawX (be sure to have the last "/" !!)
7
8 #if !defined (__CINT__) || defined (__MAKECINT__)
9
10 #include <iostream>
11 #include <stdlib.h>
12 #include <fstream>
13
14 #include "TString.h"
15 #include "TMath.h"
16
17 #include "AliHLTSystem.h"
18 #include "AliHLTPluginBase.h"
19 #include "AliLog.h"
20 #include "AliReconstruction.h"
21 #include "AliHLTDataTypes.h"
22 #include "AliHLTConfiguration.h"
23
24 #include <valgrind/callgrind.h>
25 #include <sys/time.h>
26 #include "TSystem.h"
27 #include "TFile.h"
28 #include "TGrid.h"
29 #include "TGridResult.h"
30
31 #include "AliCDBManager.h"
32 #include "AliCDBEntry.h"
33 #include "AliTriggerConfiguration.h"
34 #include "AliTriggerClass.h"
35 #include "AliExternalTrackParam.h"
36 #endif
37
38 int rec_hlt_trd(const TString input ="raw.root", TString outPath=gSystem->pwd());
39 Int_t ExtractRunNumber(const TString str);
40 int main(int argc, char** argv)
41 {
42   if(argc==2) return rec_hlt_trd(argv[1]);
43   else if(argc==3) return rec_hlt_trd(argv[1],argv[2]);
44   else return rec_hlt_trd();
45 }
46
47 int rec_hlt_trd(const TString filename, TString outPath)
48 {
49   
50   ///////////////////////////////////////////////////////////////////////////////////////////////////
51   //
52   // define the analysis chain to be run
53   //
54   
55   // What chains should be run? (usually would be: TRD-OffEsdFile)
56   TString chains="TRD-CalibFile";
57
58   // cosmics or not
59   Bool_t bCosmics=kFALSE;
60
61   // look only in data containing TRD triggers?
62   Bool_t useOnlyTRDtrigger=kFALSE;
63
64   // Is the TRD full? (this is only important if ddl files should be read)
65   Bool_t fullTRD=kTRUE;
66
67   // If not use these SMs:
68   Int_t TRDmodules[18] = {0,1,7,8,9,10,17,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1};
69
70   // Use custom arguments for components? i.e.: not reading OCDB arguments
71   Bool_t customArgs=kTRUE;
72
73   // Disable HLT flag?
74   Bool_t disableHLTflag=kFALSE;
75
76   ///////////////////////////////////////////////////////////////////////////////////////////////////
77   //
78   // init the HLT system
79   //
80   TString dataPath = outPath;
81   dataPath.Resize(dataPath.Last('/')+1);
82   if(!dataPath.Length()){
83     dataPath=gSystem->pwd();
84     dataPath+="/";
85     outPath.Prepend(dataPath);
86   }
87   
88   TString extInput = filename;
89   dataPath = filename;
90   dataPath.Resize(dataPath.Last('/')+1);
91   if(!dataPath.Length()){
92     dataPath=gSystem->pwd();
93     dataPath+="/";
94     extInput.Prepend(dataPath);
95   }
96
97   printf("File path %s\n",dataPath.Data());
98   printf("Processing file %s\n",extInput.Data());
99   printf("Output to %s\n",outPath.Data());
100
101   if(!gSystem->AccessPathName("galice.root")){
102     cerr << "please delete the galice.root or run at different place." << endl;
103     return -1;
104   }
105
106   Bool_t bRealData=kFALSE;
107   if(filename.Contains(".root") && !filename.Contains("raw.root")){
108     bRealData = kTRUE;
109     printf("processing real data\n");
110   }else{
111     bRealData = kFALSE;
112     printf("processing simulated data\n");
113   }
114
115   if(filename.Contains("alien://") || bRealData){
116     TGrid::Connect("alien://");
117   }
118
119   if(filename.Contains(".root") && !TFile::Open(filename))return -1;
120
121   gSystem->mkdir(outPath.Data());
122   gSystem->ChangeDirectory(outPath.Data());
123
124   AliHLTSystem* gHLT=AliHLTPluginBase::GetInstance();
125   
126   Int_t usedModules=0;
127   if(fullTRD){
128     usedModules = 18;
129     for(int i=0; i<18; i++)
130       TRDmodules[i]=i;
131   }else{
132     std::sort((UInt_t*)TRDmodules, ((UInt_t*)TRDmodules) + 18);
133     for(int i=0; i<18; i++)
134       if(TRDmodules[i]>-1)usedModules++;
135   }
136
137   TString option="libAliHLTUtil.so libAliHLTTRD.so libAliHLTMUON.so libAliHLTGlobal.so libAliHLTTrigger.so loglevel=0x7f chains=";
138   option+=chains;
139   TString afterTr, afterTrOff, afterCf;
140
141   for (int module = 0; module < usedModules; module++) 
142     {
143       TString arg, publisher, cf, tr, trOff;
144       // raw data publisher components
145       publisher.Form("TRD-RP_%02d", TRDmodules[module]);
146       arg.Form("-minid %d -datatype 'DDL_RAW ' 'TRD ' -dataspec %i -verbose", TRDmodules[module]+1024, (int)TMath::Power(2,TRDmodules[module]));
147       AliHLTConfiguration pubConf(publisher.Data(), "AliRawReaderPublisher", NULL , arg.Data());
148       
149       // Clusterizer
150       arg = "";
151       if(customArgs || disableHLTflag){
152         arg="output_percentage 700 -lowflux -experiment -tailcancellation -faststreamer -yPosMethod LUT";
153         if(disableHLTflag)
154           arg+=" -HLTflag no";
155       }
156
157       cf.Form("TRD-CF_%02d", TRDmodules[module]);
158       AliHLTConfiguration cfConf(cf.Data(), "TRDClusterizer", publisher.Data(), arg.Data());
159
160       if (afterCf.Length()>0) afterCf+=" ";
161       afterCf+=cf;
162
163       // Tracker
164       arg="";
165       if(customArgs || disableHLTflag){
166         arg="output_percentage 100 -lowflux -PIDmethod NN";
167         if(disableHLTflag)
168           arg+=" -HLTflag no";
169       }
170
171       tr.Form("TRD-TR_%02d", TRDmodules[module]);
172       AliHLTConfiguration trConf(tr.Data(), "TRDTrackerV1", cf.Data(), arg.Data());
173       
174       if (afterTr.Length()>0) afterTr+=" ";
175       afterTr+=tr;
176
177       // Offline Tracker (for debug purposes only)
178       arg="";
179       if(customArgs || disableHLTflag){
180         arg="output_percentage 100 -lowflux -PIDmethod NN -highLevelOutput yes -emulateHLTTracks yes";
181         if(disableHLTflag)
182           arg+=" -HLTflag no";
183       }
184
185       trOff.Form("TRD-TROFF_%02d", TRDmodules[module]);
186       AliHLTConfiguration trOffConf(trOff.Data(), "TRDOfflineTrackerV1", cf.Data(), arg.Data());
187
188       if (afterTrOff.Length()>0) afterTrOff+=" ";
189       afterTrOff+=trOff;
190       
191     }
192
193   // cluster histogramm
194   AliHLTConfiguration histoConf("TRD-ClHisto", "TRDClusterHisto", afterCf.Data(), "");
195   AliHLTConfiguration writerHistoConf( "TRD-ClHistoFile", "ROOTFileWriter", "TRD-ClHisto", "-directory hlt-trd-histo/ -datafile histo.root -concatenate-events -concatenate-blocks");
196
197   // calibration (you may use tr or trOff here)
198   AliHLTConfiguration calibConf("TRD-Calib", "TRDCalibration", afterTr.Data(), "-TrgStr hi -rejectTrgStr");
199   AliHLTConfiguration writerCalibConf( "TRD-CalibFile", "ROOTFileWriter", "TRD-Calib", "-directory hlt-trd-calib/ -datafile calib.root -concatenate-events -concatenate-blocks -write-all-events");
200
201   // esd converter
202   AliHLTConfiguration esdConf("TRD-Esd", "GlobalEsdConverter", afterTr.Data(), "-notree");
203   
204   // root file writer
205   AliHLTConfiguration writerConf("TRD-EsdFile", "EsdCollector", "TRD-Esd", "-directory hlt-trd-esd/");
206
207   // root file writer (with esd friends) (you may use tr or trOff here)
208   AliHLTConfiguration writerOffConf("TRD-OffEsdFile", "TRDEsdWriter", afterTr.Data(), "-concatenate-events -concatenate-blocks");
209
210
211   ///////////////////////////////////////////////////////////////////////////////////////////////////
212   //
213   // Init CDBManager and trigger
214   //
215   AliCDBManager * man = AliCDBManager::Instance();
216   Int_t run = ExtractRunNumber(filename);
217   if(bRealData){
218     man->SetDefaultStorage("alien://folder=/alice/data/2009/OCDB");
219   }else{
220     man->SetDefaultStorage("local://$ALICE_ROOT/OCDB"); 
221     man->SetSpecificStorage("GRP/GRP/Data", Form("local://%s",dataPath.Data()));
222     //rec.SetSpecificStorage("GRP/GRP/Data", Form("local://%s",gSystem->pwd()));
223   }
224   man->SetRun(run);
225   
226   if(bCosmics){
227     // no magnetic field
228     AliExternalTrackParam::SetMostProbablePt(8.);
229   }
230
231   // Find TRD triggers
232   TString filestring = filename;
233   if(useOnlyTRDtrigger){
234     AliCDBEntry *grp_ctp = man->Get("GRP/CTP/Config");
235     AliTriggerConfiguration *trg_conf = (AliTriggerConfiguration *)grp_ctp->GetObject();
236     trg_conf->Print();
237     TObjArray trg_masks = trg_conf->GetClasses(); // Reference!!!
238     std::vector<unsigned char> triggerconfs;
239     for(Int_t iobj = 0; iobj < trg_masks.GetEntriesFast(); iobj++){
240       
241       AliTriggerClass * trg_class = (AliTriggerClass*)trg_masks.UncheckedAt(iobj);
242       AliTriggerCluster * trg_clust = (AliTriggerCluster *)trg_class->GetCluster();
243       
244       printf("ioj[%d]\n", iobj); trg_class->Print(0x0);
245       
246       if(TString(trg_class->GetName()).Contains("TRD")){ // cosmic run 2009
247         triggerconfs.push_back(trg_class->GetMask());
248       }
249     }
250     
251     Int_t itrg = 0;
252     printf("Number of Trigger Clusters including TRD: %d\n", (Int_t)triggerconfs.size());
253     for(std::vector<unsigned char>::iterator it = triggerconfs.begin(); it < triggerconfs.end(); it++)
254       printf("Trigger Mask %d for TRD: %d\n", itrg++, *it);
255     filestring += "?EventType=7";
256     char triggerbuf[256];
257     Int_t triggerval = 0;
258     for(std::vector<unsigned char>::iterator it = triggerconfs.begin(); it < triggerconfs.end(); it++)
259       triggerval += *it;
260     sprintf(triggerbuf, "?Trigger=%d", triggerval);
261     filestring += triggerbuf; // This line does the trigger selection. It has to be uncommented if one wants to apply trigger selection
262   }
263   printf("Filename: %s\n", filestring.Data());
264
265   ///////////////////////////////////////////////////////////////////////////////////////////////////
266   //
267   // Init and run the reconstruction
268   //
269   AliReconstruction rec;
270   rec.SetInput(filestring.Data());
271   rec.SetRunVertexFinder(kFALSE);
272   rec.SetRunLocalReconstruction("HLT");
273   rec.SetRunTracking(":");
274   rec.SetLoadAlignFromCDB(0);
275   rec.SetFillESD("");
276   rec.SetRunQA(":");
277   rec.SetRunGlobalQA(kFALSE);
278   rec.SetFillTriggerESD(kFALSE);
279
280   rec.SetOption("HLT", option);
281   rec.Run();
282
283   return 0;
284 }
285
286 Int_t ExtractRunNumber(const TString str){
287   TObjArray *ptoks = (TObjArray *)str.Tokenize("?");
288   TString path = ((TObjString *)ptoks->UncheckedAt(0))->String();
289   TObjArray *toks = (TObjArray *)path.Tokenize("/");
290   TString fname = ((TObjString *)(toks->UncheckedAt(toks->GetEntriesFast() - 1)))->String();
291   TString rstr = fname(2,9);
292   return rstr.Atoi();
293 }