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