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