]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/QA/rawqa.C
typo
[u/mrichter/AliRoot.git] / test / QA / 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 "AliDAQ.h"
15 #include "AliLog.h"
16 #include "AliQA.h"
17 #include "AliQADataMakerSteer.h"
18 #include "AliRawReader.h"
19 #include "AliRawReaderRoot.h"
20 #include "AliTRDrawStreamBase.h"
21 #include "AliGeomManager.h"
22
23 TString ClassName() { return "rawqa" ; } 
24
25 //________________________________qa______________________________________
26 void rawqa(const Int_t runNumber, Int_t maxFiles = 10, const char* year = "08") 
27 {       
28         char kDefaultOCDBStorage[120] ; 
29         sprintf(kDefaultOCDBStorage, "alien://folder=/alice/data/20%s/LHC%sa/OCDB/", year, year) ; 
30         AliQA::SetQARefStorage(Form("%s%s/", AliQA::GetQARefDefaultStorage(), year)) ;  
31         AliQA::SetQARefDataDirName("Data") ; //Data, Pedestals, BlackEvent, ..... 
32         
33         UInt_t maxEvents = 99999 ;
34         if ( maxFiles < 0 ) {
35                 maxEvents = TMath::Abs(maxFiles) ; 
36                 maxFiles = 99 ; 
37         }
38         AliLog::SetGlobalDebugLevel(0) ; 
39         // connect to the grid 
40         TGrid * grid = 0x0 ; 
41         grid = TGrid::Connect("alien://") ; 
42                 
43         Bool_t detIn[AliDAQ::kNDetectors] = {kFALSE} ;
44         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"} ; 
45         // make the file name pattern year and run number
46         TString pattern;
47         pattern.Form("%9d",runNumber);
48         pattern.ReplaceAll(" ", "0") ; 
49         pattern.Prepend(year);
50         pattern.Append("*.root");
51         // find the files associated to this run
52         TGridResult * result = 0x0 ; 
53         Bool_t local = kFALSE ; 
54         if (grid) { // get the list of files from AliEn directly 
55                 TString collectionFile(pattern) ; 
56                 collectionFile.ReplaceAll("*.root", ".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 { 
61                         TString baseDir; 
62                         baseDir.Form("/alice/data/20%s/",year);
63                         result = grid->Query(baseDir, pattern) ;  
64                 }
65         } else {
66            // get the list of files from the local current directory 
67             local = kTRUE ; 
68             char line[100] ; 
69             sprintf(line, ".! ls %s*.root > tempo.txt", pattern.Data()) ; 
70             gROOT->ProcessLine(line) ;
71         }
72         
73         AliLog::Flush();
74         ifstream in ; 
75         if (local) 
76                 in.open("tempo.txt", ifstream::in) ; 
77
78         AliQADataMakerSteer qas ; 
79         TString detectors  = ""; 
80         TString detectorsW = ""; 
81         UShort_t file = 0 ; 
82         UShort_t filesProcessed = 0 ; 
83         UShort_t eventsProcessed = 0 ; 
84         AliCDBManager* man = AliCDBManager::Instance();
85         AliGeomManager::LoadGeometry("geometry.root");
86         for ( file = 0 ; file < maxFiles ; file++) {
87                 man->SetDefaultStorage(kDefaultOCDBStorage) ;  
88                 if ( qas.GetCurrentEvent() >= maxEvents) 
89                         break ;
90
91                 TString fileName ; 
92                 if ( local) {
93                         in >> fileName ;
94                 } 
95                 else 
96                         fileName = result->GetKey(file, "turl");
97                 if ( fileName == "" )  
98                         break ;
99                 if ( fileName.Contains("tag") )
100                         continue; 
101                 filesProcessed++ ;
102                 char input[200] ; 
103                 if (local) 
104                         sprintf(input, "%s", fileName.Data()) ; 
105                 else 
106                         sprintf(input, "%s", result->GetKey(file, "turl")); 
107                 AliInfo(Form("Proccessing file # %d --> %s", file, input)) ;
108                 AliLog::Flush();
109                 // check which detectors are present 
110                 AliRawReader * rawReader = new AliRawReaderRoot(input);
111                 AliTRDrawStreamBase::SetRawStreamVersion("TB");
112                 while ( rawReader->NextEvent() ) {
113                         man->SetRun(rawReader->GetRunNumber());
114                         AliLog::Flush();
115                         UChar_t * data ; 
116                         while (rawReader->ReadNextData(data)) {
117                                 Int_t detID = rawReader->GetDetectorID();
118                                 if (detID < 0 || detID >= AliDAQ::kNDetectors) {
119                                         AliError("Wrong detector ID! Skipping payload...");
120                                         continue;
121                                 }
122                                 detIn[detID] = kTRUE ; 
123                         }
124                         for (Int_t detID = 0; detID < AliDAQ::kNDetectors ; detID++) {
125                                 if (detIn[detID]) {
126                                         if ( ! detectors.Contains(detNameOff[detID]) ) {
127                                                 detectors.Append(detNameOff[detID]) ;
128                                                 detectors.Append(" ") ;
129                                         }
130                                 }
131                         }
132                         if ( !detectors.IsNull() )
133                                 break ; 
134                 }
135                 if ( !detectors.IsNull() ) {
136                         qas.SetMaxEvents(maxEvents) ;   
137                         detectorsW = qas.Run(detectors, rawReader) ;
138                         qas.Reset() ;
139                 } else {
140                         AliError("No valid detectors found") ; 
141                 } 
142                 delete rawReader ;
143                 eventsProcessed += qas.GetCurrentEvent() ; 
144         }
145         AliLog::Flush();
146         qas.Merge(runNumber) ; 
147         
148         AliLog::Flush();
149         // The summary 
150         AliInfo(Form("\n\n********** Summary for run %d **********", runNumber)) ; 
151         printf("     detectors present in the run        : %s\n", detectors.Data()) ; 
152         printf("     detectors present in the run with QA: %s\n", detectorsW.Data()) ; 
153         printf("     number of files/events processed    : %d/%d\n", filesProcessed, eventsProcessed) ; 
154         TFile * qaResult = TFile::Open(AliQA::GetQAResultFileName()) ; 
155         if ( qaResult ) {
156                 AliQA * qa = dynamic_cast<AliQA *>(qaResult->Get(AliQA::GetQAName())) ; 
157                 if ( qa) {
158                         for (Int_t index = 0 ; index < AliQA::kNDET ; index++)
159                                 if (detectorsW.Contains(AliQA::GetDetName(AliQA::DETECTORINDEX_t(index)))) 
160                                         qa->ShowStatus(AliQA::DETECTORINDEX_t(index)) ;
161                 } else {
162                         AliError(Form("%s not found in %s !", AliQA::GetQAName(), AliQA::GetQAResultFileName())) ; 
163                 }
164         } else {
165                 AliError(Form("%s has not been produced !", AliQA::GetQAResultFileName())) ; 
166         }
167 }