]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSPreprocessor.cxx
Corrected two bugs found by Sergey regarding the error calculation.
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPreprocessor.cxx
CommitLineData
1ab07e55 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
16/* $Id$ */
17
18///////////////////////////////////////////////////////////////////////////////
19// PHOS Preprocessor class. It runs by Shuttle at the end of the run,
20// calculates calibration coefficients and dead/bad channels
21// to be posted in OCDB
22//
23// Author: Boris Polichtchouk, 4 October 2006
24///////////////////////////////////////////////////////////////////////////////
25
26#include "AliPHOSPreprocessor.h"
27#include "AliLog.h"
28#include "AliCDBMetaData.h"
30a0d5a2 29#include "AliCDBEntry.h"
1ab07e55 30#include "AliPHOSEmcCalibData.h"
31#include "TFile.h"
32#include "TH1.h"
450717bc 33#include "TF1.h"
34#include "TH2.h"
1ab07e55 35#include "TMap.h"
36#include "TRandom.h"
f97eb136 37#include "TKey.h"
516e7ab3 38#include "TList.h"
39#include "TObjString.h"
30a0d5a2 40#include "TObjArray.h"
450717bc 41#include "TObject.h"
42#include "TString.h"
43#include "TMath.h"
306e5917 44#include "TAxis.h"
3c0265c1 45#include "AliPHOSEmcBadChannelsMap.h"
1ab07e55 46
47ClassImp(AliPHOSPreprocessor)
48
450717bc 49 Double_t rgaus(Double_t *x, Double_t *par)
50{
51 Double_t gaus = par[0] * TMath::Exp( -(x[0]-par[1])*(x[0]-par[1]) /
52 (2*par[2]*par[2]) );
53 return gaus;
54}
55
1ab07e55 56//_______________________________________________________________________________________
57AliPHOSPreprocessor::AliPHOSPreprocessor() :
81587b1e 58AliPreprocessor("PHS",0)
1ab07e55 59{
60 //default constructor
61}
62
63//_______________________________________________________________________________________
81587b1e 64AliPHOSPreprocessor::AliPHOSPreprocessor(AliShuttleInterface* shuttle):
65AliPreprocessor("PHS",shuttle)
1ab07e55 66{
67 // Constructor
ce2c6758 68
69 AddRunType("PHYSICS");
70 AddRunType("STANDALONE");
1ab07e55 71}
72
73//_______________________________________________________________________________________
74UInt_t AliPHOSPreprocessor::Process(TMap* /*valueSet*/)
75{
76 // process data retrieved by the Shuttle
516e7ab3 77
3c0265c1 78 TString runType = GetRunType();
79 Log(Form("Run type: %s",runType.Data()));
80
9f43d7da 81 if(runType=="STANDALONE") {
eba66a50 82 Bool_t ledOK = ProcessLEDRun();
4a7b5a26 83 Bool_t badmap_OK = FindBadChannelsEmc();
84 if(!badmap_OK) Log(Form("WARNING! FindBadChannels() completed with BAD status!"));
eba66a50 85 if(ledOK) return 0;
86 else
87 return 1;
88 }
9f43d7da 89
90 if(runType=="PHYSICS") {
91
9f43d7da 92 Bool_t calibEmc_OK = CalibrateEmc();
5db29457 93
4a7b5a26 94 if(calibEmc_OK) return 0;
9f43d7da 95 else
96 return 1;
97 }
98
99 Log(Form("Unknown run type %s. Do nothing and return OK.",runType.Data()));
100 return 0;
101
102}
103
104
105Bool_t AliPHOSPreprocessor::ProcessLEDRun()
106{
12661472 107 //Process LED run, update High Gain/Low Gain ratios.
9f43d7da 108
12661472 109 AliPHOSEmcCalibData calibData;
9f43d7da 110
111 TList* list = GetFileSources(kDAQ, "LED");
112 if(!list) {
113 Log("Sources list for LED run not found, exit.");
114 return kFALSE;
115 }
116
12661472 117 if(!list->GetEntries()) {
118 Log("Sources list for LED run is empty, exit.");
119 return kFALSE;
120 }
121
122 //Retrieve the last EMC calibration object
123 const AliPHOSEmcCalibData* clb=0;
124 AliCDBEntry* entryCalib = GetFromOCDB("Calib", "EmcGainPedestals");
125
126 if(!entryCalib)
127 Log(Form("Cannot find any AliCDBEntry for [Calib, EmcGainPedestals]!"));
128 else
129 clb = (AliPHOSEmcCalibData*)entryCalib->GetObject();
130
9f43d7da 131 TIter iter(list);
132 TObjString *source;
7708b003 133 Bool_t firedOK = kFALSE;
9f43d7da 134
135 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
12661472 136
9f43d7da 137 AliInfo(Form("found source %s", source->String().Data()));
138
139 TString fileName = GetFile(kDAQ, "LED", source->GetName());
140 AliInfo(Form("Got filename: %s",fileName.Data()));
141
142 TFile f(fileName);
143
144 if(!f.IsOpen()) {
145 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
146 return kFALSE;
147 }
148
7708b003 149 TH1I* fFiredCells = (TH1I*)f.Get("fFiredCells");
150
151 if(!fFiredCells) {
152 Log(Form("fFiredCells histogram not found in %s, skip this source.",fileName.Data()));
153 firedOK=kFALSE;
154 }
155 else {
156 const Double_t nFiredCells = fFiredCells->GetMean();
157 if(nFiredCells>100 && nFiredCells<600) {
158 firedOK = kTRUE;
159 Log(Form("Number of fired cells per event is %d.",nFiredCells));
160 }
161 else {
162 Log(Form("Number of fired cells per event is %d, possibly not a LED run!",nFiredCells));
163 firedOK = kFALSE;
164 }
165
166 }
167
9f43d7da 168 const Int_t nMod=5; // 1:5 modules
169 const Int_t nCol=56; //1:56 columns in each module
170 const Int_t nRow=64; //1:64 rows in each module
12661472 171
9f43d7da 172 for(Int_t mod=0; mod<nMod; mod++) {
173 for(Int_t col=0; col<nCol; col++) {
174 for(Int_t row=0; row<nRow; row++) {
12661472 175
7708b003 176 //High Gain to Low Gain ratio
177 if(firedOK) {
178 Float_t ratio = HG2LG(mod,row,col,&f);
179 calibData.SetHighLowRatioEmc(mod+1,col+1,row+1,ratio);
180 if(ratio != 16.)
181 AliInfo(Form("mod %d iX %d iZ %d ratio %.3f\n",mod,row,col,ratio));
182 }
183
12661472 184 if(clb) {
185 Double_t coeff = clb->GetADCchannelEmc(mod+1,col+1,row+1);
186 calibData.SetADCchannelEmc(mod+1,col+1,row+1,coeff);
187 }
188 else
189 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.);
190
9f43d7da 191 }
192 }
193 }
194
12661472 195 } // end of loop over files
9f43d7da 196
4a7b5a26 197 //Store the updated High Gain/Low Gain ratios
12661472 198 AliCDBMetaData emcMetaData;
9f43d7da 199
12661472 200 //Data valid from current run fRun until updated (validityInfinite=kTRUE)
201 Bool_t result = Store("Calib","EmcGainPedestals",&calibData,&emcMetaData,fRun,kTRUE);
9f43d7da 202 return result;
203
204}
205
206Float_t AliPHOSPreprocessor::HG2LG(Int_t mod, Int_t X, Int_t Z, TFile* f)
207{
208 //Calculates High gain to Low gain ratio
209 //for crystal at the position (X,Z) in the PHOS module mod.
210
211 char hname[128];
212 sprintf(hname,"%d_%d_%d",mod,X,Z);
213
214 TH1F* h1 = (TH1F*)f->Get(hname);
215 if(!h1) return 16.;
216
217 if(!h1->GetEntries()) return 16.;
218 if(h1->GetMaximum()<10.) return 16.;
219
220 Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
221 Double_t xmin = max - (h1->GetRMS()/3);
222 Double_t xmax = max + (h1->GetRMS()/2);
223 // Double_t xmin = max - (h1->GetRMS());
224 // Double_t xmax = max + (h1->GetRMS());
225
226 TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
227 gaus1->SetParNames("Constant","Mean","Sigma");
228 gaus1->SetParameter("Constant",h1->GetMaximum());
229 gaus1->SetParameter("Mean",max);
230 gaus1->SetParameter("Sigma",1.);
231 gaus1->SetLineColor(kBlue);
232 h1->Fit(gaus1,"LERQ+");
233
3f5f9893 234 AliInfo(Form("%s: %d entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
9f43d7da 235 h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),
236 gaus1->GetParameter("Mean")));
237
238 return gaus1->GetParameter("Mean");
239
240}
241
242Bool_t AliPHOSPreprocessor::FindBadChannelsEmc()
243{
5db29457 244 //Loop over two systems: DAQ and HLT.
245 //For each system the same algorithm implemented in DoFindBadChannelsEmc() invokes.
9f43d7da 246
5db29457 247 TList* list=0;
248 TString path;
249
250 Int_t system[2] = { kDAQ, kHLT };
251 const char* sysn[] = { "DAQ","HLT" };
252 Bool_t result[2] = { kTRUE, kTRUE };
9f43d7da 253
5db29457 254 for (Int_t i=0; i<2; i++) {
9f43d7da 255
5db29457 256 AliPHOSEmcBadChannelsMap badMap;
257 list = GetFileSources(system[i], "BAD_CHANNELS");
9f43d7da 258
5db29457 259 if(!list) {
260 Log(Form("%s sources list for BAD_CHANNELS not found!",sysn[i]));
261 result[i] = kFALSE;
262 continue;
263 }
9f43d7da 264
5db29457 265 if(!list->GetEntries()) {
266 Log(Form("Got empty sources list. It seems %s DA2 did not produce any files!",sysn[i]));
267 result[i] = kFALSE;
268 continue;
269 }
270
271 result[i] *= DoFindBadChannelsEmc(system[i],list,badMap);
272
273 // Store the bad channels map.
274
275 AliCDBMetaData md;
276 md.SetResponsible("Boris Polishchuk");
277
278 if(system[i] == kDAQ)
279 path = "Calib";
280 else
281 path = "HLT";
282
283 // Data valid from current run fRun until being updated (validityInfinite=kTRUE)
284 result[i] *= Store(path.Data(), "EmcBadChannels", &badMap, &md, fRun, kTRUE);
285
9f43d7da 286 }
5db29457 287
288 if(result[0] || result[1]) return kTRUE;
289 else return kFALSE;
290}
291
292Bool_t AliPHOSPreprocessor::DoFindBadChannelsEmc(Int_t system, TList* list, AliPHOSEmcBadChannelsMap& badMap)
293{
294 //Creates the bad channels map for PHOS EMC.
295
296 // The file fileName contains histograms which have been produced by DA2 detector algorithm.
297 // It is a responsibility of the SHUTTLE framework to form the fileName.
9f43d7da 298
299 TIter iter(list);
300 TObjString *source;
301 char hnam[80];
302 TH1F* h1=0;
303
304 const Float_t fQualityCut = 1.;
305 Int_t nGoods[5] = {0,0,0,0,0};
306
307 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
308
309 AliInfo(Form("found source %s", source->String().Data()));
310
5db29457 311 TString fileName = GetFile(system, "BAD_CHANNELS", source->GetName());
9f43d7da 312 AliInfo(Form("Got filename: %s",fileName.Data()));
313
314 TFile f(fileName);
315
316 if(!f.IsOpen()) {
317 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
318 return kFALSE;
319 }
7708b003 320
321 TH1I* fFiredCells = (TH1I*)f.Get("fFiredCells");
322
323 if(!fFiredCells) {
324 Log(Form("fFiredCells histogram not found in %s, skip this source.",fileName.Data()));
325 continue;
326 }
327
328 const Double_t nFiredCells = fFiredCells->GetMean();
329
330 if(nFiredCells<100. || nFiredCells>600. ) {
331 Log(Form("Number of fired cells per event is %d, possibly not a LED run!",nFiredCells));
332 continue; // not a LED run!
333 }
9f43d7da 334
7708b003 335 Log(Form("Number of fired cells per event %d. Begin check for bad channels.",nFiredCells));
9f43d7da 336
337 for(Int_t mod=0; mod<5; mod++) {
338 for(Int_t iX=0; iX<64; iX++) {
339 for(Int_t iZ=0; iZ<56; iZ++) {
340
341 sprintf(hnam,"%d_%d_%d_%d",mod,iX,iZ,1); // high gain
342 h1 = (TH1F*)f.Get(hnam);
343
344 if(h1) {
345 Double_t mean = h1->GetMean();
346
347 if(mean)
348 Log(Form("iX=%d iZ=%d gain=%d mean=%.3f\n",iX,iZ,1,mean));
349
350 if( mean>0 && mean<fQualityCut ) {
351 nGoods[mod]++;
352 }
353 else
354 badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
355 }
356
357 }
358 }
359
360 if(nGoods[mod])
361 Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
362 }
363
364
365 } // end of loop over sources
366
5db29457 367 return kTRUE;
9f43d7da 368}
369
370Bool_t AliPHOSPreprocessor::CalibrateEmc()
371{
5db29457 372 //Loop over two systems: DAQ and HLT.
373 //For each system the same algorithm implemented in DoCalibrateEmc() invokes.
9f43d7da 374
5db29457 375 const AliPHOSEmcBadChannelsMap* badMap=0;
376 AliCDBEntry* entryBCM=0;
377 TList* list=0;
378 TString path;
379
380 Int_t system[2] = { kDAQ, kHLT };
381 const char* sysn[] = { "DAQ","HLT" };
382 Bool_t result[2] = { kTRUE, kTRUE };
eba66a50 383
5db29457 384 for (Int_t i=0; i<2; i++) {
385
386 AliPHOSEmcCalibData calibData;
387 list = GetFileSources(system[i], "AMPLITUDES");
9f43d7da 388
5db29457 389 if(!list) {
390 Log(Form("%s sources list not found!",sysn[i]));
391 result[i] = kFALSE;
392 continue;
393 }
394
395 if(!list->GetEntries()) {
396 Log(Form("Got empty sources list. It seems %s DA1 did not produce any files!",sysn[i]));
397 result[i] = kFALSE;
398 continue;
399 }
400
401 // Retrieve the Bad Channels Map (BCM)
402
403 if(system[i] == kDAQ)
404 path = "Calib";
405 else
406 path = "HLT";
9f43d7da 407
5db29457 408 entryBCM = GetFromOCDB(path.Data(), "EmcBadChannels");
1ab07e55 409
5db29457 410 if(!entryBCM)
411 Log(Form("WARNING!! Cannot find any AliCDBEntry for [%s, EmcBadChannels]!",path.Data()));
412 else
413 badMap = (AliPHOSEmcBadChannelsMap*)entryBCM->GetObject();
414
415 if(!badMap)
416 Log(Form("WARNING!! Nothing for %s in AliCDBEntry. All cells considered GOOD!",sysn[i]));
417
418 result[i] *= DoCalibrateEmc(system[i],list,badMap,calibData);
419
420 //Store EMC calibration data
421
422 AliCDBMetaData emcMetaData;
423 result[i] *= Store(path.Data(), "EmcGainPedestals", &calibData, &emcMetaData, 0, kTRUE);
424
9f43d7da 425 }
5db29457 426
427 if(result[0] || result[1]) return kTRUE;
428 else return kFALSE;
429}
9f43d7da 430
30a0d5a2 431
5db29457 432Bool_t AliPHOSPreprocessor::DoCalibrateEmc(Int_t system, TList* list, const AliPHOSEmcBadChannelsMap* badMap, AliPHOSEmcCalibData& calibData)
433{
434 // Return kTRUE if OK.
435 // I. Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
436 // II. Extracts High_Gain/Low_Gain ratio for each channel.
30a0d5a2 437
5db29457 438 // The file fileName contains histograms which have been produced by DA1 detector algorithm.
439 // It is a responsibility of the SHUTTLE framework to form the fileName.
30a0d5a2 440
5db29457 441 gRandom->SetSeed(0); //the seed is set to the current machine clock!
30a0d5a2 442
516e7ab3 443 TIter iter(list);
444 TObjString *source;
3c0265c1 445
516e7ab3 446 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
447 AliInfo(Form("found source %s", source->String().Data()));
f97eb136 448
5db29457 449 TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
516e7ab3 450 AliInfo(Form("Got filename: %s",fileName.Data()));
f97eb136 451
516e7ab3 452 TFile f(fileName);
f97eb136 453
516e7ab3 454 if(!f.IsOpen()) {
455 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
5db29457 456 return kFALSE;
516e7ab3 457 }
9f43d7da 458
516e7ab3 459 const Int_t nMod=5; // 1:5 modules
460 const Int_t nCol=56; //1:56 columns in each module
461 const Int_t nRow=64; //1:64 rows in each module
462
463 Double_t coeff;
464 char hnam[80];
450717bc 465 TH2F* h2=0;
466 TH1D* h1=0;
3c0265c1 467
516e7ab3 468 //Get the reference histogram
469 //(author: Gustavo Conesa Balbastre)
470
471 TList * keylist = f.GetListOfKeys();
472 Int_t nkeys = f.GetNkeys();
473 Bool_t ok = kFALSE;
474 TKey *key;
516e7ab3 475 Int_t ikey = 0;
476 Int_t counter = 0;
450717bc 477 TH1D* hRef = 0;
478
516e7ab3 479 //Check if the file contains any histogram
480
481 if(nkeys< 2){
482 Log(Form("Not enough histograms (%d) for calibration.",nkeys));
4a7b5a26 483 return 1; // it's not fatal! May be short run..
f97eb136 484 }
450717bc 485
516e7ab3 486 while(!ok){
487 ikey = gRandom->Integer(nkeys);
488 key = (TKey*)keylist->At(ikey);
450717bc 489 TObject* obj = f.Get(key->GetName());
490 TString cname(obj->ClassName());
491 if(cname == "TH2F") {
492 h2 = (TH2F*)obj;
493 TString htitl = h2->GetTitle();
494 if(htitl.Contains("and gain 1")) {
495 hRef = h2->ProjectionX();
306e5917 496 hRef->GetXaxis()->SetRange(10,1000); // to cut off saturation peak and noise
450717bc 497 // Check if the reference histogram has too little statistics
306e5917 498 if(hRef->GetMean() && hRef->GetEntries()>2) ok=kTRUE;
30a0d5a2 499
500 const TString delim = "_";
501 TString str = hRef->GetName();
502 TObjArray* tks = str.Tokenize(delim);
503 const Int_t md = ((TObjString*)tks->At(0))->GetString().Atoi();
504 const Int_t X = ((TObjString*)tks->At(1))->GetString().Atoi();
505 const Int_t Z = ((TObjString*)tks->At(2))->GetString().Atoi();
506
507 if(badMap) {
508 if(badMap->IsBadChannel(md+1,Z+1,X+1)) {
509 AliInfo(Form("Cell mod=%d col=%d row=%d is bad. Histogram %s rejected.",
510 md+1,Z+1,X+1,hRef->GetName()));
511 ok=kFALSE;
512 }
513 }
514
450717bc 515 }
516e7ab3 516 }
b8093789 517
450717bc 518 counter++;
b8093789 519
520 if(!ok && counter > nkeys){
521 Log("No histogram with enough statistics for reference. Exit.");
4a7b5a26 522 return 1; // Not fatal, just wait..
b8093789 523 }
516e7ab3 524 }
3c0265c1 525
516e7ab3 526 Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
527 hRef->GetName(),hRef->GetEntries(),
528 hRef->GetMean(),hRef->GetRMS()));
529
530 Double_t refMean=hRef->GetMean();
3c0265c1 531
516e7ab3 532 // Calculates relative calibration coefficients for all non-zero channels
3c0265c1 533
516e7ab3 534 for(Int_t mod=0; mod<nMod; mod++) {
535 for(Int_t col=0; col<nCol; col++) {
536 for(Int_t row=0; row<nRow; row++) {
450717bc 537
538 sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
539 h2 = (TH2F*)f.Get(hnam);
4a7b5a26 540
516e7ab3 541 //TODO: dead channels exclusion!
450717bc 542 if(h2) {
543 h1 = h2->ProjectionX();
306e5917 544 h1->GetXaxis()->SetRange(10,1000); //to cut off saturation peak and noise
450717bc 545 coeff = h1->GetMean()/refMean;
64d8f3f1 546 if(coeff>0)
da2ee9fc 547 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1./coeff);
64d8f3f1 548 else
da2ee9fc 549 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.);
516e7ab3 550 AliInfo(Form("mod %d col %d row %d coeff %f\n",mod,col,row,coeff));
551 }
552 else
da2ee9fc 553 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.);
1ab07e55 554 }
1ab07e55 555 }
556 }
516e7ab3 557
558 f.Close();
1ab07e55 559 }
516e7ab3 560
5db29457 561 return 1;
1ab07e55 562}
450717bc 563
9f43d7da 564
450717bc 565