]> git.uio.no Git - u/mrichter/AliRoot.git/blame - test/QA/rawqa.C
typo
[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"
2dd6f761 21#include "AliGeomManager.h"
478d5b53 22
23TString ClassName() { return "rawqa" ; }
24
25//________________________________qa______________________________________
26void rawqa(const Int_t runNumber, Int_t maxFiles = 10, const char* year = "08")
27{
2dd6f761 28 AliGeomManager::LoadGeometry("geometry.root");
478d5b53 29 char kDefaultOCDBStorage[120] ;
30 sprintf(kDefaultOCDBStorage, "alien://folder=/alice/data/20%s/LHC%sa/OCDB/", year, year) ;
31 AliQA::SetQARefStorage(Form("%s%s/", AliQA::GetQARefDefaultStorage(), year)) ;
32 AliQA::SetQARefDataDirName("Data") ; //Data, Pedestals, BlackEvent, .....
33
34 UInt_t maxEvents = 99999 ;
35 if ( maxFiles < 0 ) {
36 maxEvents = TMath::Abs(maxFiles) ;
37 maxFiles = 99 ;
38 }
39 AliLog::SetGlobalDebugLevel(0) ;
40 // connect to the grid
41 TGrid * grid = 0x0 ;
42 grid = TGrid::Connect("alien://") ;
43
44 Bool_t detIn[AliDAQ::kNDetectors] = {kFALSE} ;
45 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"} ;
46 // make the file name pattern year and run number
47 TString pattern;
48 pattern.Form("%9d",runNumber);
49 pattern.ReplaceAll(" ", "0") ;
50 pattern.Prepend(year);
51 pattern.Append("*.root");
52 // find the files associated to this run
53 TGridResult * result = 0x0 ;
54 Bool_t local = kFALSE ;
55 if (grid) { // get the list of files from AliEn directly
681c8a57 56 TString collectionFile(pattern) ;
57 collectionFile.ReplaceAll("*.root", ".xml") ;
58 if ( gSystem->AccessPathName(collectionFile) == 0 ) { // get the list of files from an a-priori created collection file
59 TAlienCollection collection(collectionFile.Data(), maxFiles) ;
60 result = collection.GetGridResult("", 0, 0);
61 } else {
62 TString baseDir;
63 baseDir.Form("/alice/data/20%s/",year);
64 result = grid->Query(baseDir, pattern) ;
65 }
478d5b53 66 } else {
681c8a57 67 // get the list of files from the local current directory
478d5b53 68 local = kTRUE ;
69 char line[100] ;
70 sprintf(line, ".! ls %s*.root > tempo.txt", pattern.Data()) ;
71 gROOT->ProcessLine(line) ;
478d5b53 72 }
681c8a57 73
478d5b53 74 AliLog::Flush();
75 ifstream in ;
76 if (local)
77 in.open("tempo.txt", ifstream::in) ;
78
79 AliQADataMakerSteer qas ;
80 TString detectors = "";
81 TString detectorsW = "";
82 UShort_t file = 0 ;
83 UShort_t filesProcessed = 0 ;
84 UShort_t eventsProcessed = 0 ;
85 AliCDBManager* man = AliCDBManager::Instance();
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}