]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - PHOS/AliPHOSPreprocessor.cxx
add armenteros calculations and plots
[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 if(clb) {
164 Float_t hg2lg = clb->GetHighLowRatioEmc(5-mod,col+1,row+1);
165 Double_t coeff = clb->GetADCchannelEmc(5-mod,col+1,row+1);
166 calibData.SetADCchannelEmc(5-mod,col+1,row+1,coeff);
167 calibData.SetHighLowRatioEmc(5-mod,col+1,row+1,hg2lg);
168 }
169
170 //High Gain to Low Gain ratio
171 Float_t ratio = HG2LG(mod,row,col,&f);
172 if(ratio != 16.) {
173 calibData.SetHighLowRatioEmc(5-mod,col+1,row+1,ratio);
174 AliInfo(Form("mod %d iX %d iZ %d ratio %.3f\n",mod,row,col,ratio));
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 until updated (validityInfinite=kTRUE)
186 //Bool_t result = Store("Calib","EmcGainPedestals",&calibData,&emcMetaData,0,kTRUE);
187
188 //Store reference data
189 Bool_t refOK = StoreReferenceLED(list);
190 if(refOK) Log(Form("LED reference data successfully stored."));
191
192 //return result;
193 return kTRUE;
194}
195
196Float_t AliPHOSPreprocessor::HG2LG(Int_t mod, Int_t X, Int_t Z, TFile* f)
197{
198 //Calculates High gain to Low gain ratio
199 //for crystal at the position (X,Z) in the PHOS module mod.
200
201 char hname[128]; TString shname = "%d_%d_%d";
202 snprintf(hname,shname.Length(),shname.Data(),mod,X,Z);
203
204 TH1F* h1 = (TH1F*)f->Get(hname);
205 if(!h1) return 16.;
206
207 if(h1->GetEntries()<2000.) return 16.;
208
209 if(h1->GetMaximum()<10.) h1->Rebin(4);
210 if(h1->GetMaximum()<10.) return 16.;
211
212 Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
213 Double_t xmin = max - (h1->GetRMS()/3);
214 Double_t xmax = max + (h1->GetRMS()/2);
215 // Double_t xmin = max - (h1->GetRMS());
216 // Double_t xmax = max + (h1->GetRMS());
217
218 TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
219 gaus1->SetParNames("Constant","Mean","Sigma");
220 gaus1->SetParameter("Constant",h1->GetMaximum());
221 gaus1->SetParameter("Mean",max);
222 gaus1->SetParameter("Sigma",1.);
223 gaus1->SetLineColor(kBlue);
224
225 Double_t mean_min = h1->GetXaxis()->GetXmin();
226 Double_t mean_max = h1->GetXaxis()->GetXmax();
227 gaus1->SetParLimits(1,mean_min,mean_max);
228
229 h1->Fit(gaus1,"RQ+");
230 Double_t hg2lg = gaus1->GetParameter("Mean");
231 if( (hg2lg-mean_min<0.001) || (mean_max-hg2lg<0.001)) hg2lg=max;
232
233 AliInfo(Form("%s: %.1f entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
234 h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),hg2lg));
235
236 return hg2lg;
237
238}
239
240Bool_t AliPHOSPreprocessor::FindBadChannelsEmc()
241{
242 //Loop over two systems: DAQ and HLT.
243 //For each system the same algorithm implemented in DoFindBadChannelsEmc() invokes.
244
245 TList* list=0;
246 TString path;
247
248 Int_t system[2] = { kDAQ, kHLT };
249 const char* sysn[] = { "DAQ","HLT" };
250 Bool_t result[2] = { kTRUE, kTRUE };
251
252 for (Int_t i=0; i<2; i++) {
253
254 if(system[i] == kHLT) continue;
255
256 AliPHOSEmcBadChannelsMap badMap;
257 list = GetFileSources(system[i], "BAD_CHANNELS");
258
259 if(!list) {
260 Log(Form("%s sources list for BAD_CHANNELS not found!",sysn[i]));
261 result[i] = kFALSE;
262 continue;
263 }
264
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 Bool_t findBadOK = DoFindBadChannelsEmc(system[i],list,badMap);
272 result[i] *= findBadOK;
273
274 // Store the bad channels map.
275
276 AliCDBMetaData md;
277 md.SetResponsible("Boris Polishchuk");
278
279 if(system[i] == kDAQ)
280 path = "Calib";
281 else
282 path = "HLT";
283
284 // Data valid from current run until being updated (validityInfinite=kTRUE)
285 Bool_t storeOK = Store(path.Data(), "EmcBadChannels", &badMap, &md, 0, kTRUE);
286 result[i] *= storeOK;
287
288 }
289
290 if(result[0] || result[1]) return kTRUE;
291 else return kFALSE;
292}
293
294Bool_t AliPHOSPreprocessor::DoFindBadChannelsEmc(Int_t system, TList* list, AliPHOSEmcBadChannelsMap& badMap)
295{
296 //Creates the bad channels map for PHOS EMC.
297
298 // The file fileName contains histograms which have been produced by DA2 detector algorithm.
299 // It is a responsibility of the SHUTTLE framework to form the fileName.
300
301 TIter iter(list);
302 TObjString *source;
303 TH1F* h1=0;
304
305 const Float_t fQualityCut = 1.;
306 Int_t nGoods[5] = {0,0,0,0,0};
307
308 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
309
310 AliInfo(Form("found source %s", source->String().Data()));
311
312 TString fileName = GetFile(system, "BAD_CHANNELS", source->GetName());
313 AliInfo(Form("Got filename: %s",fileName.Data()));
314
315 TFile f(fileName);
316
317 if(!f.IsOpen()) {
318 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
319 return kFALSE;
320 }
321
322 Log(Form("Begin check for bad channels."));
323
324 for(Int_t mod=0; mod<5; mod++) {
325 for(Int_t iX=0; iX<64; iX++) {
326 for(Int_t iZ=0; iZ<56; iZ++) {
327
328 TString hnam;
329 hnam += mod; hnam += "_"; hnam += iX; hnam += "_"; hnam += iZ; hnam += "_"; hnam += "1";
330 h1 = (TH1F*)f.Get(hnam);
331
332 if(h1) {
333 Double_t mean = h1->GetMean();
334
335 if(mean)
336 Log(Form("iX=%d iZ=%d gain=%d mean=%.3f\n",iX,iZ,1,mean));
337
338 if( mean>0 && mean<fQualityCut ) {
339 nGoods[mod]++;
340 }
341 else
342 badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
343 }
344
345 }
346 }
347
348 if(nGoods[mod])
349 Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
350 }
351
352
353 } // end of loop over sources
354
355 return kTRUE;
356}
357
358Bool_t AliPHOSPreprocessor::CalibrateEmc()
359{
360 //Loop over two systems: DAQ and HLT.
361 //For each system the same algorithm implemented in DoCalibrateEmc() invokes.
362
363 AliPHOSEmcCalibData* lastCalib=0;
364 const AliPHOSEmcBadChannelsMap* badMap=0;
365 AliCDBEntry* entryBCM=0;
366 AliCDBEntry* entryEmc=0;
367 TList* list=0;
368 TString path;
369
370 Int_t system[2] = { kDAQ, kHLT };
371 const char* sysn[] = { "DAQ","HLT" };
372 Bool_t result[2] = { kTRUE, kTRUE };
373
374 for (Int_t i=0; i<2; i++) {
375
376 if(system[i] == kHLT) continue;
377
378 AliPHOSEmcCalibData calibData;
379 list = GetFileSources(system[i], "AMPLITUDES");
380
381 if(!list) {
382 Log(Form("%s sources list not found!",sysn[i]));
383 result[i] = kFALSE;
384 continue;
385 }
386
387 if(!list->GetEntries()) {
388 Log(Form("Got empty sources list. It seems %s DA1 did not produce any files!",sysn[i]));
389 result[i] = kFALSE;
390 continue;
391 }
392
393 // Retrieve the Bad Channels Map (BCM)
394
395 if(system[i] == kDAQ)
396 path = "Calib";
397 else
398 path = "HLT";
399
400 entryBCM = GetFromOCDB(path.Data(), "EmcBadChannels");
401
402 if(!entryBCM)
403 Log(Form("WARNING!! Cannot find any AliCDBEntry for [%s, EmcBadChannels]!",path.Data()));
404 else
405 badMap = (AliPHOSEmcBadChannelsMap*)entryBCM->GetObject();
406
407 if(!badMap)
408 Log(Form("WARNING!! Nothing for %s in AliCDBEntry. All cells considered GOOD!",sysn[i]));
409
410 // Retrieve the last EMC calibration object
411
412 entryEmc = GetFromOCDB(path.Data(), "EmcGainPedestals");
413
414 if(!entryEmc)
415 Log(Form("Cannot find any EmcGainPedestals entry for this run and path %s",path.Data()));
416 else
417 lastCalib = (AliPHOSEmcCalibData*)entryEmc->GetObject();
418
419 if(lastCalib)
420 result[i] *= DoCalibrateEmc(system[i],list,badMap,*lastCalib);
421 else
422 result[i] *= DoCalibrateEmc(system[i],list,badMap,calibData);
423
424 //Store EMC calibration data
425 AliCDBMetaData emcMetaData;
426
427 // if(lastCalib)
428 // result[i] *= Store(path.Data(), "EmcGainPedestals", lastCalib, &emcMetaData, 0, kFALSE);
429 // else
430 // result[i] *= Store(path.Data(), "EmcGainPedestals", &calibData, &emcMetaData, 0, kFALSE);
431
432 //Store reference data
433 Bool_t refOK = StoreReferenceEmc(system[i],list);
434 if(refOK) Log(Form("Reference data for %s amplitudes successfully stored.",sysn[i]));
435
436 }
437
438 if(result[0] || result[1]) return kTRUE;
439 else return kFALSE;
440}
441
442Bool_t AliPHOSPreprocessor::StoreReferenceEmc(Int_t system, TList* list)
443{
444 //Put 2D calibration histograms (E vs Time) prepared by DAQ/HLT to the reference storage.
445 //system is DAQ or HLT, TList is the list of FES sources.
446
447 if(system!=kDAQ) return kFALSE;
448
449 TObjString *source = dynamic_cast<TObjString *> (list->First());
450 if(!source) return kFALSE;
451
452 TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
453
454 Bool_t resultRef = StoreReferenceFile(fileName.Data(),"CalibRefPHOS.root");
455 return resultRef;
456
457}
458
459Bool_t AliPHOSPreprocessor::StoreReferenceLED(TList* list)
460{
461 //Put HG/LG histograms to the reference storage.
462
463 TObjString *source = dynamic_cast<TObjString *> (list->First());
464 if(!source) return kFALSE;
465
466 TString fileName = GetFile(kDAQ, "LED", source->GetName());
467
468 Bool_t resultRef = StoreReferenceFile(fileName.Data(),"LEDRefPHOS.root");
469 return resultRef;
470
471}
472
473
474Bool_t AliPHOSPreprocessor::DoCalibrateEmc(Int_t system, TList* list, const AliPHOSEmcBadChannelsMap* badMap, AliPHOSEmcCalibData& calibData)
475{
476 // Return kTRUE if OK.
477 // I. Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
478 // II. Extracts High_Gain/Low_Gain ratio for each channel.
479
480 // The file fileName contains histograms which have been produced by DA1 detector algorithm.
481 // It is a responsibility of the SHUTTLE framework to form the fileName.
482
483 gRandom->SetSeed(0); //the seed is set to the current machine clock!
484 Int_t minEntries=1000; // recalculate calibration coeff. if Nentries > minEntries.
485
486 TIter iter(list);
487 TObjString *source;
488
489 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
490 AliInfo(Form("found source %s", source->String().Data()));
491
492 TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
493 AliInfo(Form("Got filename: %s",fileName.Data()));
494
495 TFile f(fileName);
496
497 if(!f.IsOpen()) {
498 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
499 return kFALSE;
500 }
501
502 const Int_t nMod=5; // 1:5 modules
503 const Int_t nCol=56; //1:56 columns in each module
504 const Int_t nRow=64; //1:64 rows in each module
505
506 Double_t coeff;
507 char hnam[80];
508 TH2F* h2=0;
509 TH1D* h1=0;
510
511 //Get the reference histogram
512 //(author: Gustavo Conesa Balbastre)
513
514 TList * keylist = f.GetListOfKeys();
515 Int_t nkeys = f.GetNkeys();
516 Bool_t ok = kFALSE;
517 TKey *key;
518 Int_t ikey = 0;
519 Int_t counter = 0;
520 TH1D* hRef = 0;
521
522 //Check if the file contains any histogram
523
524 if(nkeys< 2){
525 Log(Form("Not enough histograms (%d) for calibration.",nkeys));
526 return 1; // it's not fatal! May be short run..
527 }
528
529 while(!ok){
530 ikey = gRandom->Integer(nkeys);
531 key = (TKey*)keylist->At(ikey);
532 TObject* obj = f.Get(key->GetName());
533 TString cname(obj->ClassName());
534 if(cname == "TH2F") {
535 h2 = (TH2F*)obj;
536 TString htitl = h2->GetTitle();
537 if(htitl.Contains("and gain 1")) {
538 hRef = h2->ProjectionX();
539 hRef->GetXaxis()->SetRangeUser(10.,1000.); // to cut off saturation peak and noise
540 // Check if the reference histogram has too little statistics
541 if(hRef->GetMean() && hRef->GetEntries()>minEntries) ok=kTRUE;
542
543 const TString delim = "_";
544 TString str = hRef->GetName();
545 TObjArray* tks = str.Tokenize(delim);
546 const Int_t md = ((TObjString*)tks->At(0))->GetString().Atoi();
547 const Int_t X = ((TObjString*)tks->At(1))->GetString().Atoi();
548 const Int_t Z = ((TObjString*)tks->At(2))->GetString().Atoi();
549 delete tks;
550 if(badMap) {
551 if(badMap->IsBadChannel(5-md,Z+1,X+1)) {
552 AliInfo(Form("Cell mod=%d col=%d row=%d is bad. Histogram %s rejected.",
553 5-md,Z+1,X+1,hRef->GetName()));
554 ok=kFALSE;
555 }
556 }
557
558 }
559 }
560
561 counter++;
562
563 if(!ok && counter > nkeys){
564 Log("No histogram with enough statistics for reference. Exit.");
565 return 1; // Not fatal, just wait..
566 }
567 }
568
569 Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
570 hRef->GetName(),hRef->GetEntries(),
571 hRef->GetMean(),hRef->GetRMS()));
572
573 Double_t refMean=hRef->GetMean();
574
575 // Calculates relative calibration coefficients for all non-zero channels
576
577 TString shnam = "%d_%d_%d_1";
578
579 for(Int_t mod=0; mod<nMod; mod++) {
580 for(Int_t col=0; col<nCol; col++) {
581 for(Int_t row=0; row<nRow; row++) {
582
583 snprintf(hnam,shnam.Length(),shnam.Data(),mod,row,col); // high gain!
584 h2 = (TH2F*)f.Get(hnam);
585
586 //TODO: dead channels exclusion!
587 if(h2) {
588 h1 = h2->ProjectionX();
589 h1->GetXaxis()->SetRangeUser(10.,1000.); //to cut off saturation peak and noise
590 coeff = h1->GetMean()/refMean;
591 if(coeff>0 && h1->GetEntries()>minEntries) {
592 calibData.SetADCchannelEmc(5-mod,col+1,row+1,0.005/coeff);
593 AliInfo(Form("mod %d col %d row %d coeff %f\n",mod,col,row,coeff));
594 }
595 }
596 }
597 }
598 }
599
600 f.Close();
601 }
602
603 return 1;
604}
605
606
607