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