]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/muon/TransferMuonQATrain.C
Coverity fix (an obsolete constructor removed)
[u/mrichter/AliRoot.git] / PWG3 / muon / TransferMuonQATrain.C
1 //--------------------------------------------------------------------------
2 // Macro for QA monitoring.
3 //
4 // In case it is not run with full aliroot, it needs the following libraries to compile:
5 //  - libSTEERBase.so
6 //  - libESD.so
7 //  - libAOD.so
8 //  - libANALYSIS.so
9 //  - libANALYSISalice.so
10 //  - libCORRFW.so
11 //  - libPWG3muon.so
12 //
13 //  TString includePath = "-I${ALICE_ROOT}/PWG3/base/ ";  gSystem->SetIncludePath(includePath.Data());
14
15
16 // The macro reads the PWG1 QA train output, produces a merged root files for the full period
17 // for event and track counters as well as separate root files run per run with all MUON_TRK related histograms.
18 // The results is stored under the directory "results". Then use PlotMUONQA.C, to draw QA histograms.
19 //
20 // Author: Cynthia Hadjidakis - IPN Orsay
21 //--------------------------------------------------------------------------
22
23 #if !defined(__CINT__) || defined(__MAKECINT__)
24
25 #include <Riostream.h>
26
27 // ROOT includes
28 #include "TROOT.h"
29 #include "TMath.h"
30 #include "TGrid.h"
31 #include "TGridCollection.h"
32 #include "TGridResult.h"
33 #include "TFile.h"
34 #include "TH1.h"
35 #include "TSystem.h"
36 #include "TStyle.h"
37 #include "TCanvas.h"
38 #include "TPad.h"
39
40 // ALIROOT includes
41 #include "AliCounterCollection.h"
42
43 #endif
44
45 TObjArray * GetListOfFiles(const char* baseDir, const char * trainName, const char* inFile);
46 TObjArray * GetListOfRuns(const char* runList, TObjArray *&listoffiles);
47
48 // .x TransferMuonQATrain.C("alien:///alice/data/2010/LHC10e","QA50","/Users/cynthia/Documents/alice/data/MuonQA/LHC10e/pass2/runlist_period3_test3.txt")
49 // .x TransferMuonQATrain.C("alien:///alice/cern.ch/user/s/suire/LHC10hV2116/output","","mylist.txt","AnalysisResults.root")
50 //--------------------------------------------------------------------------
51 Bool_t TransferMuonQATrain_v4(const char* baseDir, const char * trainName, const char* runList,char* inputFile = "QAresults.root")
52 {
53 #if defined(__CINT__) && !defined(__MAKECINT__)
54   gSystem->Load("libTree");
55   gSystem->Load("libGeom");
56   gSystem->Load("libVMC");
57   gSystem->Load("libPhysics");
58   gSystem->Load("libSTEERBase");
59   gSystem->Load("libESD");
60   gSystem->Load("libAOD");
61   gSystem->Load("libANALYSIS");
62   gSystem->Load("libANALYSISalice");
63   gSystem->Load("libCORRFW");
64   gSystem->Load("libPWG3base");
65   gSystem->Load("libPWG3muon");
66 #endif
67   
68   TString sbaseDir = baseDir;
69   if (sbaseDir.Contains("alien:") && !TGrid::Connect("alien://")) {
70     Error("MuonQATrain","cannot connect to grid");
71     return 0;
72   }     
73   
74   //----------------------------------------------------------- //
75   //          Build the list of files, the list of runs         //
76   //                                    to be processed
77   //----------------------------------------------------------- //
78   
79   TObjArray *listoffiles = (TObjArray*) GetListOfFiles(baseDir,trainName,inputFile);
80   if(!listoffiles) return kFALSE;
81   TObjArray *runs = (TObjArray*) GetListOfRuns(runList,listoffiles);
82   if(!runs||!listoffiles){
83     Error("TransferMuonQATrain","cannot get a list of selected runs or files");
84     return kFALSE;
85   }
86   printf("TransferMuonQATrain: Files found to be processed = %d\n",runs->GetEntriesFast());
87   if(runs->GetEntriesFast()==0) return kFALSE;
88   
89   //-------------------------------------//
90   //    Loop over the files on grid      //
91   //-------------------------------------//
92   
93   TFile *file = 0;
94   TIter next0(listoffiles);
95   TObject *nextfile;
96   TString snextfile;
97   Int_t nFiles = -1;
98   AliCounterCollection *trackCounters = 0;
99   AliCounterCollection *eventCounters = 0;
100   AliCounterCollection *mergedTrackCounters = 0;
101   AliCounterCollection *mergedEventCounters = 0;        
102   TObjArray*  outputList; 
103   TObjArray*  outputListExpert;
104   TObjArray*  outputListNorm;  
105   TFile outputHistoFile;        
106   TString fileName;
107   
108   TString command = "mkdir results";
109   cout<<"Shell command = "<<command<<endl;
110   gSystem->Exec(command.Data());
111   
112   while ((nextfile=next0())) {
113     
114     nFiles++;
115     snextfile = nextfile->GetName();
116     //Open the file
117     file = TFile::Open(snextfile.Data());
118     if(!file) continue;
119     
120     outputList = (TObjArray *)file->Get("MUON_QA/general1");
121     outputListExpert = (TObjArray *)file->Get("MUON_QA/expert");
122     outputListNorm = (TObjArray *)file->Get("MUON_QA/general2");
123     
124     trackCounters = (AliCounterCollection *) file->Get("MUON_QA/trackCounters");
125     eventCounters = (AliCounterCollection *) file->Get("MUON_QA/eventCounters");
126     
127     if(!trackCounters || !eventCounters){
128       Error("TransferMuonQATrain","Objects not found for that file");
129       continue;
130     }
131     
132     //-------------------------------------//
133     //    Merge the AliCounterCollection
134     //-------------------------------------//
135     if(nFiles==0){
136       mergedTrackCounters = (AliCounterCollection*) trackCounters->Clone();
137       mergedEventCounters = (AliCounterCollection*) eventCounters->Clone();
138     }           
139     else{
140       mergedTrackCounters->Add(trackCounters);
141       mergedEventCounters->Add(eventCounters);
142     }
143     
144     
145     //-------------------------------------//
146     //     Save the MUONQA histos in AnalysisResults.root for each run number
147     //-------------------------------------//
148     fileName = "AnalysisResults.root";
149     outputHistoFile.Open(fileName,"recreate");
150     new TDirectoryFile("MUON_QA","MUON_QA");
151     outputHistoFile.Cd("MUON_QA");
152     
153     TDirectory *dir0 = (TDirectory*) file->GetDirectory("MUON_QA");
154     TObjArray* general1 = static_cast<TObjArray*>(dir0->FindObjectAny("general1"));
155     TObjArray* expert = static_cast<TObjArray*>(dir0->FindObjectAny("expert"));
156     TObjArray* general2 = static_cast<TObjArray*>(dir0->FindObjectAny("general2"));
157     
158     if(general1) general1->Write("general1",TObject::kSingleKey);
159     if(general2) general2->Write("general2",TObject::kSingleKey);               
160     if(expert) expert->Write("expert",TObject::kSingleKey);
161     outputHistoFile.Close();
162     
163     command = "mkdir -p results/";
164     command+=((TObjString*)runs->UncheckedAt(nFiles))->GetString();
165     cout<<"Shell command = "<<command<<endl;
166     gSystem->Exec(command.Data());
167     command=" mv AnalysisResults.root results/";
168     command+=((TObjString*)runs->UncheckedAt(nFiles))->GetString();
169     command+= "/.";
170     cout<<"Shell command = "<<command<<endl;
171     gSystem->Exec(command.Data());
172   } //end of loop over files
173   
174   //-------------------------------------//
175   //      Save the AliCounterCollection in MergedAnalysisResults.root
176   //-------------------------------------//
177   
178   TFile outputFile("MergedAnalysisResults.root","recreate");
179   new TDirectoryFile("MUON_QA","MUON_QA");
180   outputFile.Cd("MUON_QA");
181   
182   mergedTrackCounters->Write();
183   mergedEventCounters->Write();
184   outputFile.Close();
185   
186   command = "mv MergedAnalysisResults.root results/.";
187   cout<<"Shell command = "<<command<<endl;
188   gSystem->Exec(command.Data());
189   
190   return kTRUE;
191 }
192
193 TObjArray * GetListOfRuns(const char* runList, TObjArray *&listoffiles)
194 {
195
196   TObjArray * runs = new TObjArray();
197   runs->SetOwner();
198   
199   if(!runList){
200     Error("GetListOfruns","runList is not defined... exit");
201     return 0;
202   }
203   else {
204     // only the ones in the runList
205     ifstream inFile(runList);
206     if (!inFile.is_open()) {
207       Error("GetListOfRuns",Form("unable to open file %s", runList));
208       return 0;
209     }
210     
211     TString currRun;
212     while (!inFile.eof()) {
213       currRun.ReadLine(inFile, kTRUE);
214       if (currRun.IsNull()) continue;
215       if (!currRun.IsDigit()) {
216         Error("GetListOfRuns","invalid run number: %s", currRun.Data());
217         return 0;
218       }
219       runs->AddLast(new TObjString(Form("%09d", currRun.Atoi())));
220     }
221     
222     inFile.close();
223   }
224         
225   printf("GetListOfRuns: Nr of runs in the runlist = %d and in the list of files = %d\n",runs->GetEntriesFast(),listoffiles->GetEntriesFast());
226   
227   if(runList && listoffiles){   
228     //Filter the selected runs and modify listoffiles
229     TObjArray*  runsFound = new TObjArray();
230     runsFound->SetOwner();      
231     
232     //filter the selected runs  
233     TIter next0(listoffiles);
234     TObject *nextfile;
235     TString snextfile;
236                 Int_t isFound = 0;
237                 
238     TObjArray *listoffilestmp = new TObjArray();        
239     listoffilestmp->SetOwner();
240     while ((nextfile=next0())) {//loop over files found on alien
241       snextfile = nextfile->GetName();
242                         isFound = 0;
243       for ( Int_t irun=0; irun<runs->GetEntriesFast(); irun++ ) { //loop over selected runs
244         TString run = ((TObjString*)runs->UncheckedAt(irun))->GetString();
245         if(snextfile.Contains(run)){
246           listoffilestmp->Add(nextfile);
247           runsFound->AddLast(new TObjString(Form("%09d", run.Atoi())));
248                 isFound=1;
249         }
250       }
251         if(isFound==0) printf("GetListOfRuns: run = %s not found in the list of files.... continue....\n",snextfile.Data());    
252
253     }
254     runs = runsFound;
255     listoffiles->Clear();
256     listoffiles = (TObjArray*) listoffilestmp->Clone();
257     
258     printf("GetListOfRuns Nr of selected runs corresponding to the list of files = %d \n",runs->GetEntriesFast());              
259   }
260   
261   return runs;
262   
263 }
264
265 TObjArray* GetListOfFiles(const char* baseDir, const char * trainName, const char* inFile)
266 {
267   
268   TString sbaseDir = baseDir;
269   TString strainName = trainName;
270   TString inputFile = inFile;
271   TString command;
272   
273   if(!sbaseDir.Contains("alien://")){
274     Error("GetListOfFiles","Not implemented for files not on alien-->exit");
275     return 0;
276   }
277   
278   sbaseDir.ReplaceAll("alien://", "");
279   
280   TObjArray *listoffiles = new TObjArray();
281   
282   if (sbaseDir.Contains(".xml")) {
283     // Read files pointed by the xml 
284     TGridCollection *coll = (TGridCollection*)gROOT->ProcessLine(Form("TAlienCollection::Open(\"%s\");", sbaseDir.Data()));
285     if (!coll) {
286       ::Error("GetListOfFiles", "Input XML collection empty.");
287       return 0;
288       }
289     // Iterate grid collection
290     while (coll->Next()) {
291       TString fname = gSystem->DirName(coll->GetTURL());
292       fname += "/";
293       fname += inputFile;      
294       listoffiles->Add(new TNamed(fname.Data(),""));
295     }   
296   }
297   else {   
298     command = Form("find %s/ *%s/%s", sbaseDir.Data(), strainName.Data(), inputFile.Data());
299     printf("command: %s\n", command.Data());
300     TGridResult *res = gGrid->Command(command);
301     if (!res) {
302       ::Error("GetListOfFiles","No result for the find command\n");
303                         cout << sbaseDir.Data() << "   " <<  strainName.Data() << endl ; 
304       delete listoffiles;
305       return 0;
306     }     
307     TIter nextmap(res);
308     TMap *map = 0;
309     while ((map=(TMap*)nextmap())) {
310       TObjString *objs = dynamic_cast<TObjString*>(map->GetValue("turl"));
311       if (!objs || !objs->GetString().Length()) {
312         // Nothing found - skip this output
313         delete res;
314         delete listoffiles;
315         return 0;
316       }
317       listoffiles->Add(new TNamed(objs->GetName(),""));
318     }
319     delete res;
320   }
321   if (!listoffiles->GetEntries()) {
322     ::Error("GetListOfFiles","No files from the find command=%s\n",command.Data());
323       delete listoffiles;
324       return 0;
325   }     
326   else printf("GetListOfFiles: Number of files found %d\n",listoffiles->GetEntries());
327   
328   return listoffiles;
329   
330 }