]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PHOS/AliPHOSPreprocessor.cxx
Added getters for ZS parameters: offset and threshold.
[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");
71eef542 70 AddRunType("LED");
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
71eef542 81 if(runType=="LED") {
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;
4d38201a 159 Log(Form("Number of fired cells per event is %.1f.",nFiredCells));
7708b003 160 }
161 else {
4d38201a 162 Log(Form("Number of fired cells per event is %.1f, possibly not a LED run!",nFiredCells));
7708b003 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
c4ed1b69 189 calibData.SetADCchannelEmc(mod+1,col+1,row+1,0.005);
12661472 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.;
1eff11b2 216
9f43d7da 217 if(!h1->GetEntries()) return 16.;
1eff11b2 218
219 h1->Rebin(4);
9f43d7da 220 if(h1->GetMaximum()<10.) return 16.;
221
222 Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
223 Double_t xmin = max - (h1->GetRMS()/3);
224 Double_t xmax = max + (h1->GetRMS()/2);
225 // Double_t xmin = max - (h1->GetRMS());
226 // Double_t xmax = max + (h1->GetRMS());
227
228 TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
229 gaus1->SetParNames("Constant","Mean","Sigma");
230 gaus1->SetParameter("Constant",h1->GetMaximum());
231 gaus1->SetParameter("Mean",max);
232 gaus1->SetParameter("Sigma",1.);
233 gaus1->SetLineColor(kBlue);
efa56161 234
235 Double_t mean_min = h1->GetXaxis()->GetXmin();
236 Double_t mean_max = h1->GetXaxis()->GetXmax();
237 gaus1->SetParLimits(1,mean_min,mean_max);
238
9f43d7da 239 h1->Fit(gaus1,"LERQ+");
240
4d38201a 241 AliInfo(Form("%s: %.1f entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
9f43d7da 242 h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),
243 gaus1->GetParameter("Mean")));
244
245 return gaus1->GetParameter("Mean");
246
247}
248
249Bool_t AliPHOSPreprocessor::FindBadChannelsEmc()
250{
5db29457 251 //Loop over two systems: DAQ and HLT.
252 //For each system the same algorithm implemented in DoFindBadChannelsEmc() invokes.
9f43d7da 253
5db29457 254 TList* list=0;
255 TString path;
256
257 Int_t system[2] = { kDAQ, kHLT };
258 const char* sysn[] = { "DAQ","HLT" };
259 Bool_t result[2] = { kTRUE, kTRUE };
9f43d7da 260
5db29457 261 for (Int_t i=0; i<2; i++) {
9f43d7da 262
5db29457 263 AliPHOSEmcBadChannelsMap badMap;
264 list = GetFileSources(system[i], "BAD_CHANNELS");
9f43d7da 265
5db29457 266 if(!list) {
267 Log(Form("%s sources list for BAD_CHANNELS not found!",sysn[i]));
268 result[i] = kFALSE;
269 continue;
270 }
9f43d7da 271
5db29457 272 if(!list->GetEntries()) {
273 Log(Form("Got empty sources list. It seems %s DA2 did not produce any files!",sysn[i]));
274 result[i] = kFALSE;
275 continue;
276 }
277
278 result[i] *= DoFindBadChannelsEmc(system[i],list,badMap);
279
280 // Store the bad channels map.
281
282 AliCDBMetaData md;
283 md.SetResponsible("Boris Polishchuk");
284
285 if(system[i] == kDAQ)
286 path = "Calib";
287 else
288 path = "HLT";
289
290 // Data valid from current run fRun until being updated (validityInfinite=kTRUE)
291 result[i] *= Store(path.Data(), "EmcBadChannels", &badMap, &md, fRun, kTRUE);
292
9f43d7da 293 }
5db29457 294
295 if(result[0] || result[1]) return kTRUE;
296 else return kFALSE;
297}
298
299Bool_t AliPHOSPreprocessor::DoFindBadChannelsEmc(Int_t system, TList* list, AliPHOSEmcBadChannelsMap& badMap)
300{
301 //Creates the bad channels map for PHOS EMC.
302
303 // The file fileName contains histograms which have been produced by DA2 detector algorithm.
304 // It is a responsibility of the SHUTTLE framework to form the fileName.
9f43d7da 305
306 TIter iter(list);
307 TObjString *source;
308 char hnam[80];
309 TH1F* h1=0;
310
311 const Float_t fQualityCut = 1.;
312 Int_t nGoods[5] = {0,0,0,0,0};
313
314 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
315
316 AliInfo(Form("found source %s", source->String().Data()));
317
5db29457 318 TString fileName = GetFile(system, "BAD_CHANNELS", source->GetName());
9f43d7da 319 AliInfo(Form("Got filename: %s",fileName.Data()));
320
321 TFile f(fileName);
322
323 if(!f.IsOpen()) {
324 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
325 return kFALSE;
326 }
7708b003 327
328 TH1I* fFiredCells = (TH1I*)f.Get("fFiredCells");
329
330 if(!fFiredCells) {
331 Log(Form("fFiredCells histogram not found in %s, skip this source.",fileName.Data()));
332 continue;
333 }
334
335 const Double_t nFiredCells = fFiredCells->GetMean();
336
337 if(nFiredCells<100. || nFiredCells>600. ) {
4d38201a 338 Log(Form("Number of fired cells per event is %.1f, possibly not a LED run!",nFiredCells));
7708b003 339 continue; // not a LED run!
340 }
9f43d7da 341
4d38201a 342 Log(Form("Number of fired cells per event %.1f. Begin check for bad channels.",nFiredCells));
9f43d7da 343
344 for(Int_t mod=0; mod<5; mod++) {
345 for(Int_t iX=0; iX<64; iX++) {
346 for(Int_t iZ=0; iZ<56; iZ++) {
347
348 sprintf(hnam,"%d_%d_%d_%d",mod,iX,iZ,1); // high gain
349 h1 = (TH1F*)f.Get(hnam);
350
351 if(h1) {
352 Double_t mean = h1->GetMean();
353
354 if(mean)
355 Log(Form("iX=%d iZ=%d gain=%d mean=%.3f\n",iX,iZ,1,mean));
356
357 if( mean>0 && mean<fQualityCut ) {
358 nGoods[mod]++;
359 }
360 else
361 badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
362 }
363
364 }
365 }
366
367 if(nGoods[mod])
368 Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
369 }
370
371
372 } // end of loop over sources
373
5db29457 374 return kTRUE;
9f43d7da 375}
376
377Bool_t AliPHOSPreprocessor::CalibrateEmc()
378{
5db29457 379 //Loop over two systems: DAQ and HLT.
380 //For each system the same algorithm implemented in DoCalibrateEmc() invokes.
9f43d7da 381
5db29457 382 const AliPHOSEmcBadChannelsMap* badMap=0;
383 AliCDBEntry* entryBCM=0;
384 TList* list=0;
385 TString path;
386
387 Int_t system[2] = { kDAQ, kHLT };
388 const char* sysn[] = { "DAQ","HLT" };
389 Bool_t result[2] = { kTRUE, kTRUE };
eba66a50 390
5db29457 391 for (Int_t i=0; i<2; i++) {
392
393 AliPHOSEmcCalibData calibData;
394 list = GetFileSources(system[i], "AMPLITUDES");
9f43d7da 395
5db29457 396 if(!list) {
397 Log(Form("%s sources list not found!",sysn[i]));
398 result[i] = kFALSE;
399 continue;
400 }
401
402 if(!list->GetEntries()) {
403 Log(Form("Got empty sources list. It seems %s DA1 did not produce any files!",sysn[i]));
404 result[i] = kFALSE;
405 continue;
406 }
407
408 // Retrieve the Bad Channels Map (BCM)
409
410 if(system[i] == kDAQ)
411 path = "Calib";
412 else
413 path = "HLT";
9f43d7da 414
5db29457 415 entryBCM = GetFromOCDB(path.Data(), "EmcBadChannels");
1ab07e55 416
5db29457 417 if(!entryBCM)
418 Log(Form("WARNING!! Cannot find any AliCDBEntry for [%s, EmcBadChannels]!",path.Data()));
419 else
420 badMap = (AliPHOSEmcBadChannelsMap*)entryBCM->GetObject();
421
422 if(!badMap)
423 Log(Form("WARNING!! Nothing for %s in AliCDBEntry. All cells considered GOOD!",sysn[i]));
424
425 result[i] *= DoCalibrateEmc(system[i],list,badMap,calibData);
426
427 //Store EMC calibration data
428
429 AliCDBMetaData emcMetaData;
430 result[i] *= Store(path.Data(), "EmcGainPedestals", &calibData, &emcMetaData, 0, kTRUE);
431
9f43d7da 432 }
5db29457 433
434 if(result[0] || result[1]) return kTRUE;
435 else return kFALSE;
436}
9f43d7da 437
30a0d5a2 438
5db29457 439Bool_t AliPHOSPreprocessor::DoCalibrateEmc(Int_t system, TList* list, const AliPHOSEmcBadChannelsMap* badMap, AliPHOSEmcCalibData& calibData)
440{
441 // Return kTRUE if OK.
442 // I. Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
443 // II. Extracts High_Gain/Low_Gain ratio for each channel.
30a0d5a2 444
5db29457 445 // The file fileName contains histograms which have been produced by DA1 detector algorithm.
446 // It is a responsibility of the SHUTTLE framework to form the fileName.
30a0d5a2 447
5db29457 448 gRandom->SetSeed(0); //the seed is set to the current machine clock!
30a0d5a2 449
516e7ab3 450 TIter iter(list);
451 TObjString *source;
3c0265c1 452
516e7ab3 453 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
454 AliInfo(Form("found source %s", source->String().Data()));
f97eb136 455
5db29457 456 TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
516e7ab3 457 AliInfo(Form("Got filename: %s",fileName.Data()));
f97eb136 458
516e7ab3 459 TFile f(fileName);
f97eb136 460
516e7ab3 461 if(!f.IsOpen()) {
462 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
5db29457 463 return kFALSE;
516e7ab3 464 }
9f43d7da 465
516e7ab3 466 const Int_t nMod=5; // 1:5 modules
467 const Int_t nCol=56; //1:56 columns in each module
468 const Int_t nRow=64; //1:64 rows in each module
469
470 Double_t coeff;
471 char hnam[80];
450717bc 472 TH2F* h2=0;
473 TH1D* h1=0;
3c0265c1 474
516e7ab3 475 //Get the reference histogram
476 //(author: Gustavo Conesa Balbastre)
477
478 TList * keylist = f.GetListOfKeys();
479 Int_t nkeys = f.GetNkeys();
480 Bool_t ok = kFALSE;
481 TKey *key;
516e7ab3 482 Int_t ikey = 0;
483 Int_t counter = 0;
450717bc 484 TH1D* hRef = 0;
485
516e7ab3 486 //Check if the file contains any histogram
487
488 if(nkeys< 2){
489 Log(Form("Not enough histograms (%d) for calibration.",nkeys));
4a7b5a26 490 return 1; // it's not fatal! May be short run..
f97eb136 491 }
450717bc 492
516e7ab3 493 while(!ok){
494 ikey = gRandom->Integer(nkeys);
495 key = (TKey*)keylist->At(ikey);
450717bc 496 TObject* obj = f.Get(key->GetName());
497 TString cname(obj->ClassName());
498 if(cname == "TH2F") {
499 h2 = (TH2F*)obj;
500 TString htitl = h2->GetTitle();
501 if(htitl.Contains("and gain 1")) {
502 hRef = h2->ProjectionX();
306e5917 503 hRef->GetXaxis()->SetRange(10,1000); // to cut off saturation peak and noise
450717bc 504 // Check if the reference histogram has too little statistics
306e5917 505 if(hRef->GetMean() && hRef->GetEntries()>2) ok=kTRUE;
30a0d5a2 506
507 const TString delim = "_";
508 TString str = hRef->GetName();
509 TObjArray* tks = str.Tokenize(delim);
510 const Int_t md = ((TObjString*)tks->At(0))->GetString().Atoi();
511 const Int_t X = ((TObjString*)tks->At(1))->GetString().Atoi();
512 const Int_t Z = ((TObjString*)tks->At(2))->GetString().Atoi();
513
514 if(badMap) {
515 if(badMap->IsBadChannel(md+1,Z+1,X+1)) {
516 AliInfo(Form("Cell mod=%d col=%d row=%d is bad. Histogram %s rejected.",
517 md+1,Z+1,X+1,hRef->GetName()));
518 ok=kFALSE;
519 }
520 }
521
450717bc 522 }
516e7ab3 523 }
b8093789 524
450717bc 525 counter++;
b8093789 526
527 if(!ok && counter > nkeys){
528 Log("No histogram with enough statistics for reference. Exit.");
4a7b5a26 529 return 1; // Not fatal, just wait..
b8093789 530 }
516e7ab3 531 }
3c0265c1 532
516e7ab3 533 Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
534 hRef->GetName(),hRef->GetEntries(),
535 hRef->GetMean(),hRef->GetRMS()));
536
537 Double_t refMean=hRef->GetMean();
3c0265c1 538
516e7ab3 539 // Calculates relative calibration coefficients for all non-zero channels
3c0265c1 540
516e7ab3 541 for(Int_t mod=0; mod<nMod; mod++) {
542 for(Int_t col=0; col<nCol; col++) {
543 for(Int_t row=0; row<nRow; row++) {
450717bc 544
545 sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
546 h2 = (TH2F*)f.Get(hnam);
4a7b5a26 547
516e7ab3 548 //TODO: dead channels exclusion!
450717bc 549 if(h2) {
550 h1 = h2->ProjectionX();
306e5917 551 h1->GetXaxis()->SetRange(10,1000); //to cut off saturation peak and noise
450717bc 552 coeff = h1->GetMean()/refMean;
64d8f3f1 553 if(coeff>0)
c4ed1b69 554 calibData.SetADCchannelEmc(mod+1,col+1,row+1,0.005/coeff);
64d8f3f1 555 else
c4ed1b69 556 calibData.SetADCchannelEmc(mod+1,col+1,row+1,0.005);
516e7ab3 557 AliInfo(Form("mod %d col %d row %d coeff %f\n",mod,col,row,coeff));
558 }
559 else
c4ed1b69 560 calibData.SetADCchannelEmc(mod+1,col+1,row+1,0.005);
1ab07e55 561 }
1ab07e55 562 }
563 }
516e7ab3 564
565 f.Close();
1ab07e55 566 }
516e7ab3 567
5db29457 568 return 1;
1ab07e55 569}
450717bc 570
9f43d7da 571
450717bc 572