]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPreprocessor.cxx
Possibility to filter the ESD friends and add objects to AliESDfriends.root. Major...
[u/mrichter/AliRoot.git] / PHOS / AliPHOSPreprocessor.cxx
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
47 ClassImp(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 //_______________________________________________________________________________________
57 AliPHOSPreprocessor::AliPHOSPreprocessor() :
58 AliPreprocessor("PHS",0)
59 {
60   //default constructor
61 }
62
63 //_______________________________________________________________________________________
64 AliPHOSPreprocessor::AliPHOSPreprocessor(AliShuttleInterface* shuttle):
65 AliPreprocessor("PHS",shuttle)
66 {
67   // Constructor
68
69   AddRunType("PHYSICS");
70   AddRunType("LED");
71 }
72
73 //_______________________________________________________________________________________
74 UInt_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
105 Bool_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           if(ratio != 16.) {
166             calibData.SetHighLowRatioEmc(5-mod,col+1,row+1,ratio);
167             AliInfo(Form("mod %d iX %d iZ %d  ratio %.3f\n",mod,row,col,ratio));
168           }
169           
170           if(clb) {
171             Double_t coeff = clb->GetADCchannelEmc(mod+1,col+1,row+1);
172             calibData.SetADCchannelEmc(mod+1,col+1,row+1,coeff);
173           }       
174         }
175       }
176     }
177
178   } // end of loop over files
179
180   //Store the updated High Gain/Low Gain ratios
181   AliCDBMetaData emcMetaData;
182
183   //Data valid from current run until updated (validityInfinite=kTRUE)
184   Bool_t result = Store("Calib","EmcGainPedestals",&calibData,&emcMetaData,0,kTRUE);
185   return result;
186
187 }
188
189 Float_t AliPHOSPreprocessor::HG2LG(Int_t mod, Int_t X, Int_t Z, TFile* f)
190 {
191   //Calculates High gain to Low gain ratio 
192   //for crystal at the position (X,Z) in the PHOS module mod.
193   
194   char hname[128];
195   sprintf(hname,"%d_%d_%d",mod,X,Z);
196
197   TH1F* h1 = (TH1F*)f->Get(hname);
198   if(!h1) return 16.;
199   
200   if(h1->GetEntries()<2000.) return 16.;
201   
202   if(h1->GetMaximum()<10.) h1->Rebin(4);
203   if(h1->GetMaximum()<10.) return 16.;
204   
205   Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
206   Double_t xmin = max - (h1->GetRMS()/3);
207   Double_t xmax = max + (h1->GetRMS()/2);
208   //       Double_t xmin = max - (h1->GetRMS());
209   //       Double_t xmax = max + (h1->GetRMS());
210
211   TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
212   gaus1->SetParNames("Constant","Mean","Sigma");
213   gaus1->SetParameter("Constant",h1->GetMaximum());
214   gaus1->SetParameter("Mean",max);
215   gaus1->SetParameter("Sigma",1.);
216   gaus1->SetLineColor(kBlue);
217   
218   Double_t mean_min = h1->GetXaxis()->GetXmin();
219   Double_t mean_max = h1->GetXaxis()->GetXmax();
220   gaus1->SetParLimits(1,mean_min,mean_max);
221   
222   h1->Fit(gaus1,"RQ+");
223   Double_t hg2lg = gaus1->GetParameter("Mean");
224   if( (hg2lg-mean_min<0.001) || (mean_max-hg2lg<0.001)) hg2lg=max;
225   
226   AliInfo(Form("%s: %.1f entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
227           h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),hg2lg)); 
228
229   return hg2lg;
230   
231 }
232
233 Bool_t AliPHOSPreprocessor::FindBadChannelsEmc()
234 {
235   //Loop over two systems: DAQ and HLT.
236   //For each system the same algorithm implemented in DoFindBadChannelsEmc() invokes.
237
238   TList* list=0;
239   TString path;
240   
241   Int_t system[2] = { kDAQ, kHLT };
242   const char* sysn[] = { "DAQ","HLT" };
243   Bool_t result[2] = { kTRUE, kTRUE };
244
245   for (Int_t i=0; i<2; i++) {
246
247     AliPHOSEmcBadChannelsMap badMap;
248     list = GetFileSources(system[i], "BAD_CHANNELS");
249
250     if(!list) {
251       Log(Form("%s sources list for BAD_CHANNELS not found!",sysn[i]));
252       result[i] = kFALSE;
253       continue;
254     }
255
256     if(!list->GetEntries()) {
257       Log(Form("Got empty sources list. It seems %s DA2 did not produce any files!",sysn[i]));
258       result[i] = kFALSE;
259       continue;
260     }
261
262     result[i] *= DoFindBadChannelsEmc(system[i],list,badMap);
263
264     // Store the bad channels map.
265   
266     AliCDBMetaData md;
267     md.SetResponsible("Boris Polishchuk");
268
269     if(system[i] == kDAQ) 
270       path = "Calib";
271     else 
272       path = "HLT";
273   
274     // Data valid from current run until being updated (validityInfinite=kTRUE)
275     result[i] *= Store(path.Data(), "EmcBadChannels", &badMap, &md, 0, kTRUE);
276     
277   }
278   
279   if(result[0] || result[1]) return kTRUE;
280   else return kFALSE;
281 }
282
283 Bool_t AliPHOSPreprocessor::DoFindBadChannelsEmc(Int_t system, TList* list, AliPHOSEmcBadChannelsMap& badMap)
284 {
285   //Creates the bad channels map for PHOS EMC.
286
287   // The file fileName contains histograms which have been produced by DA2 detector algorithm.
288   // It is a responsibility of the SHUTTLE framework to form the fileName.
289
290   TIter iter(list);
291   TObjString *source;
292   char hnam[80];
293   TH1F* h1=0;
294
295   const Float_t fQualityCut = 1.;
296   Int_t nGoods[5] = {0,0,0,0,0};
297   
298   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
299
300     AliInfo(Form("found source %s", source->String().Data()));
301
302     TString fileName = GetFile(system, "BAD_CHANNELS", source->GetName());
303     AliInfo(Form("Got filename: %s",fileName.Data()));
304
305     TFile f(fileName);
306
307     if(!f.IsOpen()) {
308       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
309       return kFALSE;
310     }
311
312     Log(Form("Begin check for bad channels."));
313
314     for(Int_t mod=0; mod<5; mod++) {
315       for(Int_t iX=0; iX<64; iX++) {
316         for(Int_t iZ=0; iZ<56; iZ++) {
317           
318           sprintf(hnam,"%d_%d_%d_%d",mod,iX,iZ,1); // high gain 
319           h1 = (TH1F*)f.Get(hnam);
320
321           if(h1) {
322             Double_t mean = h1->GetMean();
323             
324             if(mean)
325               Log(Form("iX=%d iZ=%d gain=%d   mean=%.3f\n",iX,iZ,1,mean));
326
327             if( mean>0 && mean<fQualityCut ) { 
328               nGoods[mod]++; 
329             }
330             else
331               badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
332           }
333
334         }
335       }
336
337       if(nGoods[mod])
338         Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
339     }
340
341     
342   } // end of loop over sources
343   
344   return kTRUE;
345 }
346
347 Bool_t AliPHOSPreprocessor::CalibrateEmc()
348 {
349   //Loop over two systems: DAQ and HLT.
350   //For each system the same algorithm implemented in DoCalibrateEmc() invokes.
351
352   AliPHOSEmcCalibData*   lastCalib=0;
353   const AliPHOSEmcBadChannelsMap* badMap=0;
354   AliCDBEntry* entryBCM=0;
355   AliCDBEntry* entryEmc=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     // Retrieve  the last EMC calibration object
398     
399     entryEmc = GetFromOCDB(path.Data(), "EmcGainPedestals");
400     
401     if(!entryEmc) 
402       Log(Form("Cannot find any EmcGainPedestals entry for this run and path %s",path.Data()));
403     else
404       lastCalib = (AliPHOSEmcCalibData*)entryEmc->GetObject();
405
406     if(lastCalib)
407       result[i] *= DoCalibrateEmc(system[i],list,badMap,*lastCalib);
408     else
409       result[i] *= DoCalibrateEmc(system[i],list,badMap,calibData);
410
411     //Store EMC calibration data
412
413     AliCDBMetaData emcMetaData;
414     
415     if(lastCalib)
416       result[i] *= Store(path.Data(), "EmcGainPedestals", lastCalib, &emcMetaData, 0, kTRUE);
417     else
418       result[i] *= Store(path.Data(), "EmcGainPedestals", &calibData, &emcMetaData, 0, kTRUE);
419     
420   }
421   
422   if(result[0] || result[1]) return kTRUE;
423   else return kFALSE;
424 }
425
426
427 Bool_t AliPHOSPreprocessor::DoCalibrateEmc(Int_t system, TList* list, const AliPHOSEmcBadChannelsMap* badMap, AliPHOSEmcCalibData& calibData)
428 {
429   // Return kTRUE if OK.
430   // I.  Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
431   // II. Extracts High_Gain/Low_Gain ratio for each channel.
432
433   // The file fileName contains histograms which have been produced by DA1 detector algorithm.
434   // It is a responsibility of the SHUTTLE framework to form the fileName.
435
436   gRandom->SetSeed(0); //the seed is set to the current  machine clock!
437   Int_t minEntries=100; // recalculate calibration coeff. if Nentries > minEntries.
438
439   TIter iter(list);
440   TObjString *source;
441   
442   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
443     AliInfo(Form("found source %s", source->String().Data()));
444
445     TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
446     AliInfo(Form("Got filename: %s",fileName.Data()));
447
448     TFile f(fileName);
449
450     if(!f.IsOpen()) {
451       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
452       return kFALSE;
453     }
454     
455     const Int_t nMod=5; // 1:5 modules
456     const Int_t nCol=56; //1:56 columns in each module
457     const Int_t nRow=64; //1:64 rows in each module
458
459     Double_t coeff;
460     char hnam[80];
461     TH2F* h2=0;
462     TH1D* h1=0;
463     
464     //Get the reference histogram
465     //(author: Gustavo Conesa Balbastre)
466
467     TList * keylist = f.GetListOfKeys();
468     Int_t nkeys   = f.GetNkeys();
469     Bool_t ok = kFALSE;
470     TKey  *key;
471     Int_t ikey = 0;
472     Int_t counter = 0;
473     TH1D* hRef = 0;
474
475     //Check if the file contains any histogram
476     
477     if(nkeys< 2){
478       Log(Form("Not enough histograms (%d) for calibration.",nkeys));
479       return 1; // it's not fatal! May be short run..
480     }
481     
482     while(!ok){
483       ikey = gRandom->Integer(nkeys);
484       key = (TKey*)keylist->At(ikey);
485       TObject* obj = f.Get(key->GetName());
486       TString cname(obj->ClassName());
487       if(cname == "TH2F") {
488         h2 = (TH2F*)obj;
489         TString htitl = h2->GetTitle();
490         if(htitl.Contains("and gain 1")) {
491           hRef = h2->ProjectionX();
492           hRef->GetXaxis()->SetRangeUser(10.,1000.); // to cut off saturation peak and noise
493           // Check if the reference histogram has too little statistics
494           if(hRef->GetMean() && hRef->GetEntries()>minEntries) ok=kTRUE;
495
496           const TString delim = "_";
497           TString str = hRef->GetName();
498           TObjArray* tks = str.Tokenize(delim);
499           const Int_t md = ((TObjString*)tks->At(0))->GetString().Atoi();
500           const Int_t X   = ((TObjString*)tks->At(1))->GetString().Atoi();
501           const Int_t Z   = ((TObjString*)tks->At(2))->GetString().Atoi();
502
503           if(badMap) {
504             if(badMap->IsBadChannel(5-md,Z+1,X+1)) {
505               AliInfo(Form("Cell mod=%d col=%d row=%d is bad. Histogram %s rejected.",
506                            5-md,Z+1,X+1,hRef->GetName()));
507               ok=kFALSE;
508             }
509           }
510
511         }
512       }
513       
514       counter++;
515       
516       if(!ok && counter > nkeys){
517         Log("No histogram with enough statistics for reference. Exit.");
518         return 1; // Not fatal, just wait..
519       }
520     }
521     
522     Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
523              hRef->GetName(),hRef->GetEntries(),
524              hRef->GetMean(),hRef->GetRMS()));
525
526     Double_t refMean=hRef->GetMean();
527     
528     // Calculates relative calibration coefficients for all non-zero channels
529     
530     for(Int_t mod=0; mod<nMod; mod++) {
531       for(Int_t col=0; col<nCol; col++) {
532         for(Int_t row=0; row<nRow; row++) {
533           
534           sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
535           h2 = (TH2F*)f.Get(hnam);
536           
537           //TODO: dead channels exclusion!
538           if(h2) {
539             h1 = h2->ProjectionX();
540             h1->GetXaxis()->SetRangeUser(10.,1000.); //to cut off saturation peak and noise
541             coeff = h1->GetMean()/refMean;
542             if(coeff>0 && h1->GetEntries()>minEntries) {
543               calibData.SetADCchannelEmc(5-mod,col+1,row+1,0.005/coeff);
544               AliInfo(Form("mod %d col %d row %d  coeff %f\n",mod,col,row,coeff));
545             }
546           }
547         }
548       }
549     }
550     
551     f.Close();
552   }
553   
554   return 1;
555 }
556
557      
558