]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPreprocessor.cxx
Added possibility to choose cluster time gate
[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           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
196 Float_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];
202   sprintf(hname,"%d_%d_%d",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
240 Bool_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     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 until being updated (validityInfinite=kTRUE)
284     result[i] *= Store(path.Data(), "EmcBadChannels", &badMap, &md, 0, kTRUE);
285     
286   }
287   
288   if(result[0] || result[1]) return kTRUE;
289   else return kFALSE;
290 }
291
292 Bool_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.
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
311     TString fileName = GetFile(system, "BAD_CHANNELS", source->GetName());
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     }
320
321     Log(Form("Begin check for bad channels."));
322
323     for(Int_t mod=0; mod<5; mod++) {
324       for(Int_t iX=0; iX<64; iX++) {
325         for(Int_t iZ=0; iZ<56; iZ++) {
326           
327           sprintf(hnam,"%d_%d_%d_%d",mod,iX,iZ,1); // high gain 
328           h1 = (TH1F*)f.Get(hnam);
329
330           if(h1) {
331             Double_t mean = h1->GetMean();
332             
333             if(mean)
334               Log(Form("iX=%d iZ=%d gain=%d   mean=%.3f\n",iX,iZ,1,mean));
335
336             if( mean>0 && mean<fQualityCut ) { 
337               nGoods[mod]++; 
338             }
339             else
340               badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
341           }
342
343         }
344       }
345
346       if(nGoods[mod])
347         Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
348     }
349
350     
351   } // end of loop over sources
352   
353   return kTRUE;
354 }
355
356 Bool_t AliPHOSPreprocessor::CalibrateEmc()
357 {
358   //Loop over two systems: DAQ and HLT.
359   //For each system the same algorithm implemented in DoCalibrateEmc() invokes.
360
361   AliPHOSEmcCalibData*   lastCalib=0;
362   const AliPHOSEmcBadChannelsMap* badMap=0;
363   AliCDBEntry* entryBCM=0;
364   AliCDBEntry* entryEmc=0;
365   TList* list=0;
366   TString path;
367   
368   Int_t system[2] = { kDAQ, kHLT };
369   const char* sysn[] = { "DAQ","HLT" };
370   Bool_t result[2] = { kTRUE, kTRUE };
371
372   for (Int_t i=0; i<2; i++) {
373
374     if(system[i] == kHLT) continue;
375
376     AliPHOSEmcCalibData calibData;
377     list = GetFileSources(system[i], "AMPLITUDES");
378   
379     if(!list) {
380       Log(Form("%s sources list not found!",sysn[i]));
381       result[i] = kFALSE;
382       continue;
383     }
384
385     if(!list->GetEntries()) {
386       Log(Form("Got empty sources list. It seems %s DA1 did not produce any files!",sysn[i]));
387       result[i] = kFALSE;
388       continue;
389     }
390
391     // Retrieve the Bad Channels Map (BCM)
392
393     if(system[i] == kDAQ) 
394       path = "Calib";
395     else 
396       path = "HLT";
397   
398     entryBCM = GetFromOCDB(path.Data(), "EmcBadChannels");
399
400     if(!entryBCM)
401       Log(Form("WARNING!! Cannot find any AliCDBEntry for [%s, EmcBadChannels]!",path.Data()));
402     else
403       badMap = (AliPHOSEmcBadChannelsMap*)entryBCM->GetObject();
404
405     if(!badMap)
406       Log(Form("WARNING!! Nothing for %s in AliCDBEntry. All cells considered GOOD!",sysn[i]));
407
408     // Retrieve  the last EMC calibration object
409     
410     entryEmc = GetFromOCDB(path.Data(), "EmcGainPedestals");
411     
412     if(!entryEmc) 
413       Log(Form("Cannot find any EmcGainPedestals entry for this run and path %s",path.Data()));
414     else
415       lastCalib = (AliPHOSEmcCalibData*)entryEmc->GetObject();
416
417     if(lastCalib) 
418       result[i] *= DoCalibrateEmc(system[i],list,badMap,*lastCalib);    
419     else 
420       result[i] *= DoCalibrateEmc(system[i],list,badMap,calibData);
421     
422     //Store EMC calibration data
423     AliCDBMetaData emcMetaData;
424     
425     // if(lastCalib)
426     //   result[i] *= Store(path.Data(), "EmcGainPedestals", lastCalib, &emcMetaData, 0, kFALSE);
427     // else
428     //   result[i] *= Store(path.Data(), "EmcGainPedestals", &calibData, &emcMetaData, 0, kFALSE);
429
430     //Store reference data
431     Bool_t refOK = StoreReferenceEmc(system[i],list);
432     if(refOK) Log(Form("Reference data for %s amplitudes successfully stored.",sysn[i]));
433     
434   }
435   
436   if(result[0] || result[1]) return kTRUE;
437   else return kFALSE;
438 }
439
440 Bool_t AliPHOSPreprocessor::StoreReferenceEmc(Int_t system, TList* list)
441 {
442   //Put 2D calibration histograms (E vs Time) prepared by DAQ/HLT to the reference storage.
443   //system is DAQ or HLT, TList is the list of FES sources.
444
445   if(system!=kDAQ) return kFALSE;
446
447   TObjString *source = dynamic_cast<TObjString *> (list->First());
448   if(!source) return kFALSE;
449
450   TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
451
452   Bool_t resultRef = StoreReferenceFile(fileName.Data(),"CalibRefPHOS.root");
453   return resultRef;
454
455 }
456
457 Bool_t AliPHOSPreprocessor::StoreReferenceLED(TList* list)
458 {
459   //Put HG/LG histograms to the reference storage.
460   
461   TObjString *source = dynamic_cast<TObjString *> (list->First());
462   if(!source) return kFALSE;
463   
464   TString fileName = GetFile(kDAQ, "LED", source->GetName());
465   
466   Bool_t resultRef = StoreReferenceFile(fileName.Data(),"LEDRefPHOS.root");
467   return resultRef;
468   
469 }
470
471
472 Bool_t AliPHOSPreprocessor::DoCalibrateEmc(Int_t system, TList* list, const AliPHOSEmcBadChannelsMap* badMap, AliPHOSEmcCalibData& calibData)
473 {
474   // Return kTRUE if OK.
475   // I.  Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
476   // II. Extracts High_Gain/Low_Gain ratio for each channel.
477
478   // The file fileName contains histograms which have been produced by DA1 detector algorithm.
479   // It is a responsibility of the SHUTTLE framework to form the fileName.
480
481   gRandom->SetSeed(0); //the seed is set to the current  machine clock!
482   Int_t minEntries=1000; // recalculate calibration coeff. if Nentries > minEntries.
483
484   TIter iter(list);
485   TObjString *source;
486   
487   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
488     AliInfo(Form("found source %s", source->String().Data()));
489
490     TString fileName = GetFile(system, "AMPLITUDES", source->GetName());
491     AliInfo(Form("Got filename: %s",fileName.Data()));
492
493     TFile f(fileName);
494
495     if(!f.IsOpen()) {
496       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
497       return kFALSE;
498     }
499     
500     const Int_t nMod=5; // 1:5 modules
501     const Int_t nCol=56; //1:56 columns in each module
502     const Int_t nRow=64; //1:64 rows in each module
503
504     Double_t coeff;
505     char hnam[80];
506     TH2F* h2=0;
507     TH1D* h1=0;
508     
509     //Get the reference histogram
510     //(author: Gustavo Conesa Balbastre)
511
512     TList * keylist = f.GetListOfKeys();
513     Int_t nkeys   = f.GetNkeys();
514     Bool_t ok = kFALSE;
515     TKey  *key;
516     Int_t ikey = 0;
517     Int_t counter = 0;
518     TH1D* hRef = 0;
519
520     //Check if the file contains any histogram
521     
522     if(nkeys< 2){
523       Log(Form("Not enough histograms (%d) for calibration.",nkeys));
524       return 1; // it's not fatal! May be short run..
525     }
526     
527     while(!ok){
528       ikey = gRandom->Integer(nkeys);
529       key = (TKey*)keylist->At(ikey);
530       TObject* obj = f.Get(key->GetName());
531       TString cname(obj->ClassName());
532       if(cname == "TH2F") {
533         h2 = (TH2F*)obj;
534         TString htitl = h2->GetTitle();
535         if(htitl.Contains("and gain 1")) {
536           hRef = h2->ProjectionX();
537           hRef->GetXaxis()->SetRangeUser(10.,1000.); // to cut off saturation peak and noise
538           // Check if the reference histogram has too little statistics
539           if(hRef->GetMean() && hRef->GetEntries()>minEntries) ok=kTRUE;
540
541           const TString delim = "_";
542           TString str = hRef->GetName();
543           TObjArray* tks = str.Tokenize(delim);
544           const Int_t md = ((TObjString*)tks->At(0))->GetString().Atoi();
545           const Int_t X   = ((TObjString*)tks->At(1))->GetString().Atoi();
546           const Int_t Z   = ((TObjString*)tks->At(2))->GetString().Atoi();
547
548           if(badMap) {
549             if(badMap->IsBadChannel(5-md,Z+1,X+1)) {
550               AliInfo(Form("Cell mod=%d col=%d row=%d is bad. Histogram %s rejected.",
551                            5-md,Z+1,X+1,hRef->GetName()));
552               ok=kFALSE;
553             }
554           }
555
556         }
557       }
558       
559       counter++;
560       
561       if(!ok && counter > nkeys){
562         Log("No histogram with enough statistics for reference. Exit.");
563         return 1; // Not fatal, just wait..
564       }
565     }
566     
567     Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
568              hRef->GetName(),hRef->GetEntries(),
569              hRef->GetMean(),hRef->GetRMS()));
570
571     Double_t refMean=hRef->GetMean();
572     
573     // Calculates relative calibration coefficients for all non-zero channels
574     
575     for(Int_t mod=0; mod<nMod; mod++) {
576       for(Int_t col=0; col<nCol; col++) {
577         for(Int_t row=0; row<nRow; row++) {
578           
579           sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
580           h2 = (TH2F*)f.Get(hnam);
581           
582           //TODO: dead channels exclusion!
583           if(h2) {
584             h1 = h2->ProjectionX();
585             h1->GetXaxis()->SetRangeUser(10.,1000.); //to cut off saturation peak and noise
586             coeff = h1->GetMean()/refMean;
587             if(coeff>0 && h1->GetEntries()>minEntries) {
588               calibData.SetADCchannelEmc(5-mod,col+1,row+1,0.005/coeff);
589               AliInfo(Form("mod %d col %d row %d  coeff %f\n",mod,col,row,coeff));
590             }
591           }
592         }
593       }
594     }
595     
596     f.Close();
597   }
598   
599   return 1;
600 }
601
602      
603