2 contact: Boris.Polishchuk@cern.ch
\r
3 reference run: /alice/data/2010/LHC10a_PHOS/000112189/raw/10000112189003.10.root
\r
6 number of events needed: 200
\r
7 input files: Mod0RCU0.data Mod0RCU1.data Mod0RCU2.data Mod0RCU3.data Mod1RCU0.data Mod1RCU1.data Mod1RCU2.data Mod1RCU3.data Mod2RCU0.data Mod2RCU1.data Mod2RCU2.data Mod2RCU3.data Mod3RCU0.data Mod3RCU1.data Mod3RCU2.data Mod3RCU3.data Mod4RCU0.data Mod4RCU1.data Mod4RCU2.data Mod4RCU3.data
\r
8 Output files: PHOS_PED.root
\r
13 #include "monitor.h"
\r
22 #include <TSystem.h>
\r
24 #include <TPluginManager.h>
\r
30 #include <TString.h>
\r
31 #include "AliRawReader.h"
\r
32 #include "AliCaloRawStreamV3.h"
\r
35 int main(int argc, char **argv)
\r
37 // Read raw data, decode it to samples,
\r
38 // calculate pedestals from presamples,
\r
39 // evaluate the signal amplitude as a maximum sample,
\r
40 // and fill histograms with pedestals and amplitudes
\r
42 gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
\r
48 AliLog::SetGlobalDebugLevel(0) ;
\r
49 AliLog::SetGlobalLogLevel(AliLog::kFatal);
\r
54 printf("Wrong number of arguments\n");
\r
58 /* Retrieve mapping files from DAQ DB */
\r
59 const char* mapFiles[20] = {
\r
82 for(Int_t iFile=0; iFile<20; iFile++) {
\r
83 int failed = daqDA_DB_getFile(mapFiles[iFile], mapFiles[iFile]);
\r
85 printf("Cannot retrieve file %s from DAQ DB. Exit.\n",mapFiles[iFile]);
\r
90 /* Open mapping files */
\r
91 AliAltroMapping *mapping[20];
\r
92 TString path = "./";
\r
99 for(Int_t iMod = 0; iMod < 5; iMod++) {
\r
104 for(Int_t iRCU=0; iRCU<4; iRCU++) {
\r
108 mapping[iMap] = new AliCaloAltroMapping(path3.Data());
\r
113 /* define data source : this is argument 1 */
\r
114 status=monitorSetDataSource( argv[1] );
\r
116 printf("monitorSetDataSource() failed : %s\n",monitorDecodeError(status));
\r
120 /* declare monitoring program */
\r
121 status=monitorDeclareMp( __FILE__ );
\r
123 printf("monitorDeclareMp() failed : %s\n",monitorDecodeError(status));
\r
127 /* define wait event timeout - 1s max */
\r
128 monitorSetNowait();
\r
129 monitorSetNoWaitNetworkTimeout(1000);
\r
131 /* init some counters */
\r
132 int nevents_physics=0;
\r
133 int nevents_total=0;
\r
135 AliRawReader * reader = NULL;
\r
136 AliCaloRawStreamV3* stream = NULL;
\r
138 TString baseNamePed ="hPed";
\r
139 TString baseTitlePed="Ped in cell (";
\r
140 const char* sgain[3]={"LG","HG", "TRU"};
\r
142 const Int_t caloFlagMax=3,modMax=5,cellXMax=64,cellZMax=56;
\r
143 TH1F *hPed[5][3][64][56] = {};
\r
145 TH1F *hPedHiMean1m2 = new TH1F("hPedHiMean1m2","Mean pedestals in module 2, high gain" ,100,0.,100.);
\r
146 TH1F *hPedHiRMS1m2 = new TH1F("hPedHiRMS1m2" ,"RMS pedestals in module 2, high gain" ,100,0.,50.);
\r
147 TH1F *hPedLoMean1m2 = new TH1F("hPedLoMean1m2","Mean pedestals in module 2, low gain" ,100,0.,100.);
\r
148 TH1F *hPedLoRMS1m2 = new TH1F("hPedLoRMS1m2" ,"RMS pedestals in module 2, low gain" ,100,0.,50.);
\r
149 TH1F *hPedTRUMean1m2 = new TH1F("hPedTRUMean1m2","Mean pedestals in module 2, TRU" ,1000,0.,1000.);
\r
150 TH1F *hPedTRURMS1m2 = new TH1F("hPedTRURMS1m2" ,"RMS pedestals in module 2, TRU" ,100,0.,50.);
\r
152 TH1F *hPedHiMean1m3 = new TH1F("hPedHiMean1m3","Mean pedestals in module 3, high gain" ,100,0.,100.);
\r
153 TH1F *hPedHiRMS1m3 = new TH1F("hPedHiRMS1m3" ,"RMS pedestals in module 3, high gain" ,100,0.,50.);
\r
154 TH1F *hPedLoMean1m3 = new TH1F("hPedLoMean1m3","Mean pedestals in module 3, low gain" ,100,0.,100.);
\r
155 TH1F *hPedLoRMS1m3 = new TH1F("hPedLoRMS1m3" ,"RMS pedestals in module 3, low gain" ,100,0.,50.);
\r
156 TH1F *hPedTRUMean1m3 = new TH1F("hPedTRUMean1m3","Mean pedestals in module 3, TRU" ,1000,0.,1000.);
\r
157 TH1F *hPedTRURMS1m3 = new TH1F("hPedTRURMS1m3" ,"RMS pedestals in module 3, TRU" ,100,0.,50.);
\r
159 TH1F *hPedHiMean1m4 = new TH1F("hPedHiMean1m4","Mean pedestals in module 4, high gain" ,100,0.,100.);
\r
160 TH1F *hPedHiRMS1m4 = new TH1F("hPedHiRMS1m4" ,"RMS pedestals in module 4, high gain" ,100,0.,50.);
\r
161 TH1F *hPedLoMean1m4 = new TH1F("hPedLoMean1m4","Mean pedestals in module 4, low gain" ,100,0.,100.);
\r
162 TH1F *hPedLoRMS1m4 = new TH1F("hPedLoRMS1m4" ,"RMS pedestals in module 4, low gain" ,100,0.,50.);
\r
163 TH1F *hPedTRUMean1m4 = new TH1F("hPedTRUMean1m4","Mean pedestals in module 4, TRU" ,1000,0.,1000.);
\r
164 TH1F *hPedTRURMS1m4 = new TH1F("hPedTRURMS1m4" ,"RMS pedestals in module 4, TRU" ,100,0.,50.);
\r
166 hPedHiMean1m2->Sumw2();
\r
167 hPedHiRMS1m2 ->Sumw2();
\r
168 hPedLoMean1m2->Sumw2();
\r
169 hPedLoRMS1m2 ->Sumw2();
\r
170 hPedTRUMean1m2->Sumw2();
\r
171 hPedTRURMS1m2 ->Sumw2();
\r
172 hPedHiMean1m3->Sumw2();
\r
173 hPedHiRMS1m3 ->Sumw2();
\r
174 hPedLoMean1m3->Sumw2();
\r
175 hPedLoRMS1m3 ->Sumw2();
\r
176 hPedTRUMean1m3->Sumw2();
\r
177 hPedTRURMS1m3 ->Sumw2();
\r
178 hPedHiMean1m4->Sumw2();
\r
179 hPedHiRMS1m4 ->Sumw2();
\r
180 hPedLoMean1m4->Sumw2();
\r
181 hPedLoRMS1m4 ->Sumw2();
\r
182 hPedTRUMean1m4->Sumw2();
\r
183 hPedTRURMS1m4 ->Sumw2();
\r
185 TH2F *hPedHiMeanm2 = new TH2F("hPedHiMeanm2","Mean pedestals in module 2, high gain",
\r
186 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
187 TH2F *hPedHiRMSm2 = new TH2F("hPedHiRMSm2" ,"R.M.S. of pedestals in module 2, high gain",
\r
188 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
189 TH2F *hPedHiNumm2 = new TH2F("hPedHiNumm2" ,"Number of pedestals in module 2, high gain",
\r
190 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
191 TH2F *hPedLoMeanm2 = new TH2F("hPedLoMeanm2","Mean pedestals in module 2, low gain",
\r
192 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
193 TH2F *hPedLoRMSm2 = new TH2F("hPedLoRMSm2" ,"R.M.S. of pedestals in module 2, low gain",
\r
194 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
195 TH2F *hPedLoNumm2 = new TH2F("hPedLoNumm2" ,"Number of pedestals in module 2, low gain",
\r
196 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
198 TH2F *hPedHiMeanm3 = new TH2F("hPedHiMeanm3","Mean pedestals in module 3, high gain",
\r
199 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
200 TH2F *hPedHiRMSm3 = new TH2F("hPedHiRMSm3" ,"R.M.S. of pedestals in module 3, high gain",
\r
201 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
202 TH2F *hPedHiNumm3 = new TH2F("hPedHiNumm3" ,"Number of pedestals in module 3, high gain",
\r
203 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
204 TH2F *hPedLoMeanm3 = new TH2F("hPedLoMeanm3","Mean pedestals in module 3, low gain",
\r
205 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
206 TH2F *hPedLoRMSm3 = new TH2F("hPedLoRMSm3" ,"R.M.S. of pedestals in module 3, low gain",
\r
207 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
208 TH2F *hPedLoNumm3 = new TH2F("hPedLoNumm3" ,"Number of pedestals in module 3, low gain",
\r
209 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
211 TH2F *hPedHiMeanm4 = new TH2F("hPedHiMeanm4","Mean pedestals in module 4, high gain",
\r
212 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
213 TH2F *hPedHiRMSm4 = new TH2F("hPedHiRMSm4" ,"R.M.S. of pedestals in module 4, high gain",
\r
214 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
215 TH2F *hPedHiNumm4 = new TH2F("hPedHiNumm4" ,"Number of pedestals in module 4, high gain",
\r
216 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
217 TH2F *hPedLoMeanm4 = new TH2F("hPedLoMeanm4","Mean pedestals in module 4, low gain",
\r
218 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
219 TH2F *hPedLoRMSm4 = new TH2F("hPedLoRMSm4" ,"R.M.S. of pedestals in module 4, low gain",
\r
220 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
221 TH2F *hPedLoNumm4 = new TH2F("hPedLoNumm4" ,"Number of pedestals in module 4, low gain",
\r
222 cellXMax,0.,cellXMax, cellZMax,0.,cellZMax);
\r
224 TH1I *hNBunches = new TH1I("hNBunches","Number of bunches",10,0,10);
\r
225 TH2I *hHWaddr = new TH2I("hHWaddr","DDL is vs HW addr",216,0,216,4096,0,4096);
\r
226 TH1I *hModule = new TH1I("hModule" ,"Module number", 5,0.,5);
\r
229 Int_t module,cellX,cellZ,caloFlag;
\r
231 /* main loop (infinite) */
\r
233 struct eventHeaderStruct *event;
\r
234 eventTypeType eventT;
\r
236 /* check shutdown condition */
\r
237 if (daqDA_checkShutdown()) {break;}
\r
239 /* get next event (blocking call until timeout) */
\r
240 status=monitorGetEventDynamic((void **)&event);
\r
241 if (status==MON_ERR_EOF) {
\r
242 printf ("End of File detected\n");
\r
243 break; /* end of monitoring file has been reached */
\r
247 printf("monitorGetEventDynamic() failed : %s\n",monitorDecodeError(status));
\r
251 /* retry if got no event */
\r
256 /* use event - here, just write event id to result file */
\r
257 eventT=event->eventType;
\r
259 if (eventT==PHYSICS_EVENT) {
\r
261 reader = new AliRawReaderDate((void*)event);
\r
262 stream = new AliCaloRawStreamV3(reader,"PHOS",mapping);
\r
264 while (reader->NextEvent()) {
\r
265 runNum = reader->GetRunNumber();
\r
266 while (stream->NextDDL()) {
\r
267 while (stream->NextChannel()) {
\r
268 module = stream->GetModule();
\r
269 cellX = stream->GetCellX();
\r
270 cellZ = stream->GetCellZ();
\r
271 caloFlag = stream->GetCaloFlag();
\r
272 if (caloFlag!=0 && caloFlag!=1) continue;
\r
274 hHWaddr->Fill(stream->GetDDLNumber(),stream->GetHWAddress());
\r
275 hModule->Fill(module);
\r
276 if (!hPed[module][caloFlag][cellX][cellZ]) {
\r
277 TString name = baseNamePed;
\r
278 TString title = baseTitlePed;
\r
279 name +="_g"; name +=caloFlag;
\r
280 name +="_m"; name +=module;
\r
281 name +="_x"; name +=cellX;
\r
282 name +="_z"; name +=cellZ;
\r
284 title +=module; title +=",";
\r
285 title +=cellX; title +=",";
\r
286 title +=cellZ; title +="), ";
\r
287 title +=sgain[caloFlag];
\r
289 Int_t nx,xmin,xmax;
\r
290 if (caloFlag==0 || caloFlag==1) {
\r
300 hPed[module][caloFlag][cellX][cellZ] = new TH1F(name,title,100,0.,100.);
\r
301 hPed[module][caloFlag][cellX][cellZ]->Sumw2();
\r
302 hPed[module][caloFlag][cellX][cellZ]->SetMarkerStyle(20);
\r
303 hPed[module][caloFlag][cellX][cellZ]->SetOption("eph");
\r
306 Int_t nBunches = 0;
\r
307 while (stream->NextBunch()) {
\r
309 const UShort_t *sig = stream->GetSignals();
\r
310 Int_t sigLength = stream->GetBunchLength();
\r
311 for (Int_t i = 0; i < sigLength; i++) {
\r
312 hPed[module][caloFlag][cellX][cellZ]->Fill(sig[i]);
\r
315 hNBunches->Fill(nBunches);
\r
316 } // end of NextChannel()
\r
318 } // end of NextDDL()
\r
319 } // end of nextEvent()
\r
321 // Fill 2-dim histograms for mean, rms and n pedestals
\r
323 for (Int_t mod=2; mod<=4; mod++) {
\r
324 for (Int_t caloFlag=0; caloFlag<2; caloFlag++) {
\r
325 for (Int_t cellX=0; cellX<cellXMax; cellX++) {
\r
326 for (Int_t cellZ=0; cellZ<cellZMax; cellZ++) {
\r
327 if (hPed[mod][caloFlag][cellX][cellZ] != 0) {
\r
328 if (caloFlag == 0) {
\r
330 hPedLoMean1m2->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
331 hPedLoRMS1m2 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
332 hPedLoMeanm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
333 hPedLoRMSm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
334 hPedLoNumm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());
\r
337 hPedLoMean1m3->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
338 hPedLoRMS1m3 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
339 hPedLoMeanm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
340 hPedLoRMSm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
341 hPedLoNumm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());
\r
344 hPedLoMean1m4->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
345 hPedLoRMS1m4 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
346 hPedLoMeanm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
347 hPedLoRMSm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
348 hPedLoNumm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());
\r
351 else if (caloFlag == 1) {
\r
353 hPedHiMean1m2->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
354 hPedHiRMS1m2 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );
\r
355 hPedHiMeanm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
356 hPedHiRMSm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
357 hPedHiNumm2 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());
\r
360 hPedHiMean1m3->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
361 hPedHiRMS1m3 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );
\r
362 hPedHiMeanm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
363 hPedHiRMSm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
364 hPedHiNumm3 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());
\r
367 hPedHiMean1m4->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
368 hPedHiRMS1m4 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
369 hPedHiMeanm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
370 hPedHiRMSm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetRMS());
\r
371 hPedHiNumm4 ->Fill( cellX, cellZ, hPed[mod][caloFlag][cellX][cellZ]->GetEntries());
\r
374 else if (caloFlag == 2) {
\r
376 hPedTRUMean1m2->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
377 hPedTRURMS1m2 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );
\r
380 hPedTRUMean1m3->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
381 hPedTRURMS1m3 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );
\r
384 hPedTRUMean1m4->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetMean());
\r
385 hPedTRURMS1m4 ->Fill( hPed[mod][caloFlag][cellX][cellZ]->GetRMS() );
\r
397 } // end of if (eventT==PHYSICS_EVENT)
\r
401 /* free resources */
\r
404 /* exit when last event received, no need to wait for TERM signal */
\r
405 if (eventT==END_OF_RUN) {
\r
406 printf("EOR event detected\n");
\r
410 } // end of inf. loop over events
\r
412 for(Int_t i = 0; i < 20; i++) delete mapping[i];
\r
414 // Write existing histograms to a root file
\r
416 TString fileName = "PHOS_PED.root";
\r
417 TFile *file = new TFile(fileName,"RECREATE");
\r
419 for (Int_t mod=2; mod<=3; mod++) {
\r
420 for (Int_t caloFlag=0; caloFlag<caloFlagMax; caloFlag++) {
\r
421 for (Int_t mod=0; mod<modMax; mod++) {
\r
422 for (Int_t cellX=0; cellX<cellXMax; cellX++) {
\r
423 for (Int_t cellZ=0; cellZ<cellZMax; cellZ++) {
\r
424 if (hPed[mod][caloFlag][cellX][cellZ] != 0)
\r
425 hPed[mod][caloFlag][cellX][cellZ]->Write();
\r
432 hPedHiMean1m2->Write();
\r
433 hPedHiRMS1m2 ->Write();
\r
434 hPedLoMean1m2->Write();
\r
435 hPedLoRMS1m2 ->Write();
\r
436 hPedHiMeanm2 ->Write();
\r
437 hPedHiRMSm2 ->Write();
\r
438 hPedHiNumm2 ->Write();
\r
439 hPedLoMeanm2 ->Write();
\r
440 hPedLoRMSm2 ->Write();
\r
441 hPedLoNumm2 ->Write();
\r
442 hPedTRUMean1m2->Write();
\r
443 hPedTRURMS1m2 ->Write();
\r
445 hPedHiMean1m3->Write();
\r
446 hPedHiRMS1m3 ->Write();
\r
447 hPedLoMean1m3->Write();
\r
448 hPedLoRMS1m3 ->Write();
\r
449 hPedHiMeanm3 ->Write();
\r
450 hPedHiRMSm3 ->Write();
\r
451 hPedHiNumm3 ->Write();
\r
452 hPedLoMeanm3 ->Write();
\r
453 hPedLoRMSm3 ->Write();
\r
454 hPedLoNumm3 ->Write();
\r
455 hPedTRUMean1m3->Write();
\r
456 hPedTRURMS1m3 ->Write();
\r
458 hPedHiMean1m4->Write();
\r
459 hPedHiRMS1m4 ->Write();
\r
460 hPedLoMean1m4->Write();
\r
461 hPedLoRMS1m4 ->Write();
\r
462 hPedHiMeanm4 ->Write();
\r
463 hPedHiRMSm4 ->Write();
\r
464 hPedHiNumm4 ->Write();
\r
465 hPedLoMeanm4 ->Write();
\r
466 hPedLoRMSm4 ->Write();
\r
467 hPedLoNumm4 ->Write();
\r
468 hPedTRUMean1m4->Write();
\r
469 hPedTRURMS1m4 ->Write();
\r
471 hNBunches ->Write();
\r
477 /* Store output files to the File Exchange Server */
\r
478 daqDA_FES_storeFile(fileName.Data(),"PED");
\r
480 printf("%d physics events of %d total processed.\n",nevents_physics,nevents_total);
\r