]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/macros/makeResults.C
7c46c9c5c9eb365c948f7f85dd4947edd3f4ca80
[u/mrichter/AliRoot.git] / TRD / qaRec / macros / makeResults.C
1 // Usage:
2 //   makeResults.C("tasks", "file_list", kGRID)
3 //   tasks : "ALL" or one/more of the following separated by space:
4 //     "EFF"  : TRD Tracking Efficiency 
5 //     "EFFC" : TRD Tracking Efficiency Combined (barrel + stand alone) - only in case of simulations
6 //     "RES"  : TRD tracking Resolution
7 //     "CAL"  : TRD calibration
8 //     "PID"  : TRD PID - pion efficiency 
9 //     "PIDR" : TRD PID - reference data
10 //     "DET"  : Basic TRD Detector checks
11 //     "NOFR" : Data set does not have AliESDfriends.root 
12 //     "NOMC" : Data set does not have Monte Carlo Informations (real data), so all tasks which rely
13 //              on MC information are switched off
14 //   file_list : is the list of the files to be processed. 
15 //     They should contain the full path. Here is an example:
16 // /lustre/alice/local/TRDdata/SIM/P-Flat/TRUNK/RUN0/TRD.Performance.root
17 // or for GRID alien:///alice/cern.ch/user/m/mfasel/MinBiasProd/results/ppMinBias80000/1/TRD.Performance.root
18 //   kGRID : specify if files are comming from a GRID collection [default kFALSE]
19 //
20 // HOW TO BUILD THE FILE LIST
21 //   1. locally
22 // ls -1 BaseDir/RUN*/TRD.Performance.root > files.lst
23 // 
24 //   2. on Grid
25 // char *BaseDir="/alice/cern.ch/user/m/mfasel/MinBiasProd/results/ppMinBias80000/";
26 // TGrid::Connect("alien://");
27 // TGridResult *res = gGrid->Query(BaseDir, "%/TRD.Task%.root");
28 // TGridCollection *col = gGrid->OpenCollectionQuery(res);
29 // col->Reset();
30 // TMap *map = 0x0;
31 // while(map = (TMap*)col->Next()){
32 //   TIter it((TCollection*)map);
33 //   TObjString *info = 0x0;
34 //   while(info=(TObjString*)it()){
35 //     printf("alien://%s\n", col->GetLFN(info->GetString().Data()));
36 //   }
37 // }
38 //
39 // The files which will be processed are the intersection between the
40 // condition on the tasks and the files in the file list.
41 //
42 // Authors:
43 //   Alex Bercuci (A.Bercuci@gsi.de) 
44 //   Markus Fasel (m.Fasel@gsi.de) 
45 //
46
47 #if ! defined (__CINT__) || defined (__MAKECINT__)
48 #include <fstream>
49 #include "TError.h"
50 #include <TClass.h>
51 #include <TFileMerger.h>
52 #include <TCanvas.h>
53 #include <TH1.h>
54 #include <TGraph.h>
55 #include <TObjArray.h>
56 #include <TObjString.h>
57 #include <TString.h>
58 #include <TROOT.h>
59 #include <TSystem.h>
60 #include <TStyle.h>
61
62 #include "qaRec/AliTRDrecoTask.h"
63
64 #endif
65
66 #include "macros/AliTRDperformanceTrain.h"
67 //#include "../../PWG1/macros/AddPerformanceTask.h"
68
69 Char_t *libs[] = {"libProofPlayer.so", "libANALYSIS.so", "libTRDqaRec.so"};
70
71 void mergeProd(const Char_t *mark="TRD.Performance.root", const Char_t *files=0);
72 void makeResults(Char_t *opt = "ALL", const Char_t *files=0x0, Bool_t kGRID=kFALSE)
73 {
74   if(kGRID){
75     if(!gSystem->Getenv("GSHELL_ROOT")){
76       Error("makeResults.C", "AliEn not initialized.");
77       return;
78     }
79     TGrid::Connect("alien://");
80   }
81
82         // Load Libraries in interactive mode
83   Int_t nlibs = static_cast<Int_t>(sizeof(libs)/sizeof(Char_t *));
84   for(Int_t ilib=0; ilib<nlibs; ilib++){
85     if(gSystem->Load(libs[ilib]) >= 0) continue;
86     Error("makeResults.C", Form("Failed to load %s.", libs[ilib]));
87     return;
88   }
89
90   // define setup 
91   gStyle->SetOptStat(0);
92   Bool_t mc      = kTRUE;
93   Bool_t friends = kTRUE;
94
95
96   if(files) mergeProd("TRD.Performance.root", files);
97   Int_t fSteerTask = ParseOptions(opt);
98   TCanvas *c=new TCanvas("c", "TRD Performance", 10, 10, 800, 500);
99
100   TClass *ctask = new TClass;
101   AliTRDrecoTask *task = 0x0;
102   for(Int_t itask = NTRDQATASKS; itask--;){
103     if(!TSTBIT(fSteerTask, itask)) continue;
104
105     new(ctask) TClass(fgkTRDtaskClassName[itask]);
106     task = (AliTRDrecoTask*)ctask->New();
107     task->SetDebugLevel(0);
108     task->SetMCdata(mc);
109     task->SetFriends(friends);
110
111     if(!task->Load(Form("%s/TRD.Performance.root", gSystem->ExpandPathName("$PWD")))){
112       Error("makeResults.C", Form("Load data container for task %s failed.", task->GetName()));
113       delete task;
114       break;
115     }
116
117     if(!task->PostProcess()){
118       Error("makeResults.C", Form("Processing data container for task %s failed.", task->GetName()));
119       delete task;
120       break;
121     }
122     for(Int_t ipic=0; ipic<task->GetNRefFigures(); ipic++){
123       c->Clear();
124       if(!task->GetRefFigure(ipic)) continue;
125       c->SaveAs(Form("%s_Fig%02d.gif", task->GetName(), ipic));
126     }
127     delete task;
128   }
129   delete ctask;
130   delete c;
131 }
132
133 //______________________________________________________
134 void mergeProd(const Char_t *mark, const Char_t *files)
135 {
136   const Int_t kBatch = 20;
137
138   TFileMerger *fFM = new TFileMerger(1);
139   fFM->OutputFile(Form("%s/0_%s",  gSystem->ExpandPathName("$PWD"), mark));
140
141   Int_t jbatch = 0, nbatch = 0;
142   string filename;
143   ifstream file(files);
144   while(getline(file, filename)){
145     if(Int_t(filename.find(mark)) < 0) continue;
146     fFM->AddFile(filename.c_str()); nbatch++;
147     if(nbatch==kBatch){
148       //printf("MERGE BATCH %d [%d]\n", jbatch, nbatch);
149       fFM->Merge(); jbatch++;
150       new(fFM) TFileMerger(kTRUE);
151       fFM->OutputFile(Form("%s/%d_%s",  gSystem->ExpandPathName("$PWD"), jbatch, mark));
152       nbatch=0;
153     }
154   }
155   if(nbatch){
156     //printf("MERGE BATCH %d[%d]\n", jbatch, nbatch);
157     fFM->Merge();
158     jbatch++;
159   }
160   if(!jbatch){
161     delete fFM;
162     return;
163   }
164
165   new(fFM) TFileMerger(kTRUE);
166   fFM->OutputFile(Form("%s/%s",  gSystem->ExpandPathName("$PWD"), mark));
167   for(Int_t ib=jbatch; ib--;){
168     fFM->AddFile(Form("%s/%d_%s",  gSystem->ExpandPathName("$PWD"), ib, mark));
169     gSystem->Exec(Form("rm -f %s/%d_%s", gSystem->ExpandPathName("$PWD"), ib, mark));
170   }
171   fFM->Merge();
172   delete fFM;
173 }
174