]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/QA/rawqa.C
TuneA settings removed from kPyJets
[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 = Form("alien://folder=/alice/data/20%s/LHC%sd/OCDB/", year, year) ; 
29         char * kDefaultOCDBStorage = Form("local://$ALICE_ROOT/OCDB") ; 
30         //AliQA::SetQARefStorage(Form("%s%s/", AliQA::GetQARefDefaultStorage(), year)) ;  
31         AliQA::SetQARefStorage("local://$ALICE_ROOT/QAref") ;
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         TString 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         TString detectors  = ""; 
79         TString detectorsW = ""; 
80         UShort_t file = 0 ; 
81         UShort_t filesProcessed = 0 ; 
82         UShort_t eventsProcessed = 0 ; 
83         AliCDBManager* man = AliCDBManager::Instance();
84         man->SetDefaultStorage(kDefaultOCDBStorage) ;  
85         man->SetRun(runNumber) ; 
86         AliQAManager qam("rec") ; 
87   qam->SetDefaultStorage(AliQA::GetQARefStorage()) ; 
88   qam->SetRun(runNumber) ; 
89     
90         AliGeomManager::LoadGeometry();
91         for ( file = 0 ; file < maxFiles ; file++) {
92                 if ( qas.GetCurrentEvent() >= maxEvents) 
93                         break ;
94
95                 TString fileName ; 
96                 if ( local) {
97                         in >> fileName ;
98                 } 
99                 else 
100                         fileName = result->GetKey(file, "turl");
101                 if ( fileName == "" )  
102                         break ;
103                 if ( fileName.Contains("tag") )
104                         continue; 
105                 filesProcessed++ ;
106                 char input[200] ; 
107                 if (local) 
108                         sprintf(input, "%s", fileName.Data()) ; 
109                 else 
110                         sprintf(input, "%s", result->GetKey(file, "turl")); 
111                 AliInfo(Form("Proccessing file # %d --> %s", file, input)) ;
112                 AliLog::Flush();
113                 // check which detectors are present 
114                 AliRawReader * rawReader = new AliRawReaderRoot(input);
115                 AliTRDrawStreamBase::SetRawStreamVersion("TB");
116                 while ( rawReader->NextEvent() ) {
117                         man->SetRun(rawReader->GetRunNumber());
118                         qam->SetRun(rawReader->GetRunNumber());
119                         AliLog::Flush();
120                         UChar_t * data ; 
121                         while (rawReader->ReadNextData(data)) {
122                                 Int_t detID = rawReader->GetDetectorID();
123                                 if (detID < 0 || detID >= AliDAQ::kNDetectors) {
124                                         AliError("Wrong detector ID! Skipping payload...");
125                                         continue;
126                                 }
127                                 detIn[detID] = kTRUE ; 
128                         }
129                         for (Int_t detID = 0; detID < AliDAQ::kNDetectors ; detID++) {
130                                 if (detIn[detID]) {
131                                         if ( ! detectors.Contains(detNameOff[detID]) ) {
132                                                 detectors.Append(detNameOff[detID]) ;
133                                                 detectors.Append(" ") ;
134                                         }
135                                 }
136                         }
137                         if ( !detectors.IsNull() )
138                                 break ; 
139                 }
140                 if ( !detectors.IsNull() ) {
141                         qam.SetMaxEvents(maxEvents) ;   
142                         detectorsW = qam.Run(detectors, rawReader) ;
143                         qam.Reset() ;
144                 } else {
145                         AliError("No valid detectors found") ; 
146                 } 
147                 delete rawReader ;
148                 eventsProcessed += qas.GetCurrentEvent() ; 
149         }
150         AliLog::Flush();
151         qam.Merge(runNumber) ; 
152         
153         AliLog::Flush();
154         // The summary 
155         AliInfo(Form("\n\n********** Summary for run %d **********", runNumber)) ; 
156         printf("     detectors present in the run        : %s\n", detectors.Data()) ; 
157         printf("     detectors present in the run with QA: %s\n", detectorsW.Data()) ; 
158         printf("     number of files/events processed    : %d/%d\n", filesProcessed, eventsProcessed) ; 
159         TFile * qaResult = TFile::Open(AliQA::GetQAResultFileName()) ; 
160         if ( qaResult ) {
161                 AliQA * qa = dynamic_cast<AliQA *>(qaResult->Get(AliQA::GetQAName())) ; 
162                 if ( qa) {
163                         for (Int_t index = 0 ; index < AliQA::kNDET ; index++)
164                                 if (detectorsW.Contains(AliQA::GetDetName(AliQA::DETECTORINDEX_t(index)))) 
165                                         qa->ShowStatus(AliQA::DETECTORINDEX_t(index)) ;
166                 } else {
167                         AliError(Form("%s not found in %s !", AliQA::GetQAName(), AliQA::GetQAResultFileName())) ; 
168                 }
169         } else {
170                 AliError(Form("%s has not been produced !", AliQA::GetQAResultFileName())) ; 
171         }
172 }