]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/cosmic/rawqa.C
added max events to process
[u/mrichter/AliRoot.git] / test / cosmic / rawqa.C
1 #include <iostream>
2 #include <fstream>
3
4 #include <TAlienCollection.h>
5 #include <TFile.h>
6 #include <TGrid.h>
7 #include <TGridResult.h>
8 #include <TMath.h>
9 #include <TROOT.h>
10 #include <TString.h>
11 #include <TSystem.h>
12
13 #include "AliCDBManager.h"
14 #include "AliLog.h"
15 #include "AliQA.h"
16 #include "AliQADataMakerSteer.h"
17 #include "AliRawReader.h"
18 #include "AliRawReaderRoot.h"
19 #include "AliDAQ.h"
20
21 TString ClassName() { return "rawqa" ; } 
22
23 //________________________________qa______________________________________
24 void rawqa(const Int_t runNumber, Int_t maxFiles = 10, const char* year = "08") 
25 {       
26         const char * kDefaultOCDBStorage = Form("alien://folder=/alice/data/20%s/LHC%sa/OCDB/", year, year) ; 
27         
28         UInt_t maxEvents = 99999 ;
29         if ( maxFiles < 0 ) {
30                 maxEvents = TMath::Abs(maxFiles) ; 
31                 maxFiles = 99 ; 
32         }
33         AliLog::SetGlobalDebugLevel(0) ; 
34         // connect to the grid 
35         TGrid * grid = 0x0 ; 
36         grid = TGrid::Connect("alien://") ; 
37                 
38         Bool_t detIn[AliDAQ::kNDetectors] = {kFALSE} ;
39         char * detNameOff[AliDAQ::kNDetectors] = {"ITS", "ITS", "ITS", "TPC", "TRD", "TOF", "HMPID", "PHOS", "PHOS", "PMD", "MUON", "MUON", "FMD", "T0", "VZERO", "ZDC", "ACORDE", "TRG", "EMCAL", "DAQ_TEST", "HLT"} ; 
40         // make the file name pattern year and run number
41         TString pattern;
42         pattern.Form("%9d",runNumber);
43         pattern.ReplaceAll(" ", "0") ; 
44         pattern.Prepend(year);
45         pattern.Append("*0.root");
46
47         // find the files associated to this run
48         TGridResult * result = 0x0 ; 
49         Bool_t local = kFALSE ; 
50         if (grid) { // get the list of files from AliEn directly 
51           TString baseDir; 
52           baseDir.Form("/alice/data/20%s/",year);
53           result = grid->Query(baseDir, pattern) ;  
54         } else {
55           TString collectionFile(pattern) ; 
56           collectionFile.Append(".xml") ; 
57           if ( gSystem->AccessPathName(collectionFile) == 0 ) { // get the list of files from an a-priori created collection file
58             TAlienCollection collection(collectionFile.Data(), maxFiles) ; 
59             result = collection.GetGridResult("", 0, 0); 
60           } else { // get the list of files from the local current directory 
61             local = kTRUE ; 
62             char line[100] ; 
63             sprintf(line, ".! ls %s*.root > tempo.txt", pattern.Data()) ; 
64             gROOT->ProcessLine(line) ;
65           }
66         }
67         AliLog::Flush();
68         ifstream in ; 
69         if (local) 
70                 in.open("tempo.txt", ifstream::in) ; 
71
72         AliCDBManager* man = AliCDBManager::Instance();
73         man->SetDefaultStorage(kDefaultOCDBStorage) ;  
74         AliQA::SetQARefStorage("alien://folder=/alice/QA/2008/") ; 
75         man->SetSpecificStorage("*", AliQA::GetQARefStorage());
76         AliQADataMakerSteer qas ; 
77         TString detectors  = ""; 
78         TString detectorsW = ""; 
79         UShort_t file = 0 ; 
80         UShort_t filesProcessed = 0 ; 
81         for ( file = 0 ; file < maxFiles ; file++) {
82                 TString fileName ; 
83                 if ( local) {
84                         in >> fileName ;
85                 } 
86                 else 
87                         fileName = result->GetKey(file, "turl");
88                 if ( fileName == "" )  
89                         break ;
90                 if ( fileName.Contains("tag") )
91                         continue; 
92                 filesProcessed++ ;
93                 char input[200] ; 
94                 if (local) 
95                         sprintf(input, "%s", fileName.Data()) ; 
96                 else 
97                         sprintf(input, "%s", result->GetKey(file, "turl")); 
98                 AliInfo(Form("Proccessing file # %d --> %s", file, input)) ;
99                 AliLog::Flush();
100                 // check which detectors are present 
101                 AliRawReader * rawReader = new AliRawReaderRoot(input);
102                 while ( rawReader->NextEvent() ) {
103                         man->SetRun(rawReader->GetRunNumber());
104                         AliLog::Flush();
105                         UChar_t * data ; 
106                         while (rawReader->ReadNextData(data)) {
107                                 Int_t detID = rawReader->GetDetectorID();
108                                 if (detID < 0 || detID >= AliDAQ::kNDetectors) {
109                                         AliError("Wrong detector ID! Skipping payload...");
110                                         continue;
111                                 }
112                                 detIn[detID] = kTRUE ; 
113                         }
114                         for (Int_t detID = 0; detID < AliDAQ::kNDetectors ; detID++) {
115                                 if (detIn[detID]) {
116                                         if ( ! detectors.Contains(detNameOff[detID]) ) {
117                                                 detectors.Append(detNameOff[detID]) ;
118                                                 detectors.Append(" ") ;
119                                         }
120                                 }
121                         }
122                         if ( !detectors.IsNull() )
123                                 break ; 
124                 }
125                 AliInfo(Form("Current = %d Max = %d", qas.GetCurrentEvent(), maxEvents)) ; 
126                 AliLog::Flush();
127                 if ( qas.GetCurrentEvent() > maxEvents) 
128                         break ;
129                 // TEMPORARY REMOVAL OF TRD!!!
130                 detectors.ReplaceAll("TRD", "") ;
131                 // TEMPORARY REMOVAL OF TRD!!!
132                 if ( !detectors.IsNull() ) {
133                         qas.SetMaxEvents(maxEvents) ;   
134                         detectorsW = qas.Run(detectors, rawReader) ;
135                         qas.Reset() ;
136                 } else {
137                         AliError("No valid detectors found") ; 
138                 } 
139                 delete rawReader ;
140         }
141         AliLog::Flush();
142         qas.Merge(runNumber) ; 
143         
144         AliLog::Flush();
145         // The summary 
146         AliInfo(Form("\n\n********** Summary for run %d **********", runNumber)) ; 
147         printf("     detectors present in the run        : %s\n", detectors.Data()) ; 
148         printf("     detectors present in the run with QA: %s\n", detectorsW.Data()) ; 
149         printf("     number of files processed           : %d\n", filesProcessed) ; 
150         TFile * qaResult = TFile::Open(AliQA::GetQAResultFileName()) ; 
151         AliQA * qa = dynamic_cast<AliQA *>(qaResult->Get("QA")) ; 
152         for (Int_t index = 0 ; index < AliQA::kNDET ; index++)
153                 if (detectorsW.Contains(AliQA::GetDetName(AliQA::DETECTORINDEX(index)))) 
154                         qa->ShowStatus(AliQA::DETECTORINDEX(index)) ;
155 }