]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/macros/rawmerge.C
- Macro to create the "raw" data file with selected events
[u/mrichter/AliRoot.git] / PWG1 / macros / rawmerge.C
1 //\r
2 // Macro to create the "raw" data file with selected events\r
3 // Paramaters:\r
4 //  eventListFileName - input ascii file with list of chunks and event numbers - two collumns\r
5 //  osplit            - split the file after given threhshold value\r
6 //\r
7 void rawmerge(const char *eventListFileName,\r
8               const char *outputDirectoryURI,\r
9               Long64_t osplit=-1)\r
10 {\r
11   Int_t eventNumber;\r
12   FILE *files=fopen(eventListFileName,"r");\r
13   if (!files) {\r
14     fprintf(stderr,"error: could not read event list file \"%s\". Exiting.\n",eventListFileName);\r
15     return;\r
16   }\r
17   char iURI[1000];\r
18   char oURI[1000];\r
19   TFile *ifile=0;\r
20   TFile *ofile=0;\r
21   TTree *itree=0;\r
22   TTree *otree=0;\r
23   Long64_t ievent;\r
24   Long64_t oevent;\r
25   Int_t ofilenumber=0;\r
26   Int_t line=0;\r
27   while (!feof(files)) {\r
28     delete ifile;ifile=0;\r
29     ++line;\r
30     if (fscanf(files,"%s %d\n",iURI,&ievent)!=2) {\r
31       fprintf(stderr,"warning: corrupted event line (%d) in input file, skipping it...\n",line);\r
32       continue;\r
33     }\r
34     printf("> processing \"%s\" event %d...\n",iURI,ievent);\r
35     TFile *ifile=TFile::Open(iURI);\r
36     if (!ifile) {\r
37       fprintf(stderr,"warning: could not open file for event \"%s\", skipping it...\n",iURI);\r
38       continue;\r
39     }\r
40     TTree *itree=dynamic_cast<TTree*>(ifile->Get("RAW"));\r
41     if (!itree) {\r
42       fprintf(stderr,"warning: could not find RAW tree for event \"%s\", skipping it...\n",iURI);\r
43       continue;\r
44     }\r
45 \r
46     // create (new) output file and tree\r
47     if (!ofile || (osplit>0 && oevent%osplit==0)) {\r
48       delete ofile;\r
49       ++ofilenumber;\r
50       sprintf(oURI,"%s/merged_%d.root",outputDirectoryURI,ofilenumber);\r
51       printf("< creating output file \"%s\"\n",oURI);\r
52       ofile=TFile::Open(oURI,"RECREATE");\r
53       if (!ofile) {\r
54         fprintf(stderr,"error: could not create output file: \"%s\" Exiting.\n",oURI);\r
55         break;\r
56       }\r
57       otree=itree->CloneTree(0);\r
58     }\r
59 \r
60     // copy event and write to file\r
61     otree->CopyAddresses(itree);\r
62     itree->GetEntry(ievent);\r
63     otree->Fill();\r
64 //    otree->CopyEntries(itree,Form("Entry$==%d",ievent),1);\r
65     ofile->Write();\r
66     ++oevent;\r
67 \r
68     // reset input\r
69     itree->ResetBranchAddresses();\r
70   }\r
71 \r
72   printf("Merged %d events.\n",oevent);\r
73   delete ifile;\r
74   delete ofile;\r
75 }\r
76 \r