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