]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSPreprocessor.cxx
Histogram ranges changed to cut off saturation peak and noise
[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 "AliPHOSEmcCalibData.h"
30 #include "TFile.h"
31 #include "TH1.h"
32 #include "TF1.h"
33 #include "TH2.h"
34 #include "TMap.h"
35 #include "TRandom.h"
36 #include "TKey.h"
37 #include "TList.h"
38 #include "TObjString.h"
39 #include "TObject.h"
40 #include "TString.h"
41 #include "TMath.h"
42 #include "TAxis.h"
43 #include "AliPHOSEmcBadChannelsMap.h"
44
45 ClassImp(AliPHOSPreprocessor)
46
47   Double_t rgaus(Double_t *x, Double_t *par)
48 {
49   Double_t gaus = par[0] * TMath::Exp( -(x[0]-par[1])*(x[0]-par[1]) /
50                                        (2*par[2]*par[2]) );
51   return gaus;
52 }
53
54 //_______________________________________________________________________________________
55 AliPHOSPreprocessor::AliPHOSPreprocessor() :
56 AliPreprocessor("PHS",0)
57 {
58   //default constructor
59 }
60
61 //_______________________________________________________________________________________
62 AliPHOSPreprocessor::AliPHOSPreprocessor(AliShuttleInterface* shuttle):
63 AliPreprocessor("PHS",shuttle)
64 {
65   // Constructor
66 }
67
68 //_______________________________________________________________________________________
69 UInt_t AliPHOSPreprocessor::Process(TMap* /*valueSet*/)
70 {
71   // process data retrieved by the Shuttle
72   
73   TString runType = GetRunType();
74   Log(Form("Run type: %s",runType.Data()));
75
76   if(runType=="STANDALONE") {
77     Bool_t ledOK = ProcessLEDRun();
78     if(ledOK) return 0;
79     else
80       return 1;
81   }
82   
83   if(runType=="PHYSICS") {
84     
85     Bool_t badmap_OK = FindBadChannelsEmc();
86     if(!badmap_OK) Log(Form("WARNING!! FindBadChannels() completed with BAD status!"));
87     
88     Bool_t calibEmc_OK = CalibrateEmc();
89
90     if(calibEmc_OK && badmap_OK) return 0;
91     else
92       return 1;
93   }
94
95   Log(Form("Unknown run type %s. Do nothing and return OK.",runType.Data()));
96   return 0;
97
98 }
99
100
101 Bool_t AliPHOSPreprocessor::ProcessLEDRun()
102 {
103   //Process LED run, fill bad channels map.
104
105   AliPHOSEmcBadChannelsMap badMap;
106
107   TList* list = GetFileSources(kDAQ, "LED");
108   if(!list) {
109     Log("Sources list for LED run not found, exit.");
110     return kFALSE;
111   }
112
113   TIter iter(list);
114   TObjString *source;
115   char hnam[80];
116   TH1F* histo=0;
117   
118   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
119
120     AliInfo(Form("found source %s", source->String().Data()));
121
122     TString fileName = GetFile(kDAQ, "LED", source->GetName());
123     AliInfo(Form("Got filename: %s",fileName.Data()));
124
125     TFile f(fileName);
126
127     if(!f.IsOpen()) {
128       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
129       return kFALSE;
130     }
131
132     const Int_t nMod=5; // 1:5 modules
133     const Int_t nCol=56; //1:56 columns in each module
134     const Int_t nRow=64; //1:64 rows in each module
135
136     // Check for dead channels    
137     Log(Form("Begin check for dead channels."));
138
139     for(Int_t mod=0; mod<nMod; mod++) {
140       for(Int_t col=0; col<nCol; col++) {
141         for(Int_t row=0; row<nRow; row++) {
142           sprintf(hnam,"mod%dcol%drow%d",mod,col,row);
143           histo = (TH1F*)f.Get(hnam);
144           if(histo)
145             if (histo->GetMean()<1) {
146               Log(Form("Channel: [%d,%d,%d] seems dead, <E>=%.1f.",mod,col,row,histo->GetMean()));
147               badMap.SetBadChannel(mod,col,row);
148             }
149         }
150       }
151     }
152
153   }
154
155   //Store bad channels map
156   AliCDBMetaData badMapMetaData;
157
158   //Bad channels data valid from current run fRun until updated (validityInfinite=kTRUE)
159   Bool_t result = Store("Calib", "EmcBadChannels", &badMap, &badMapMetaData, fRun, kTRUE);
160   return result;
161
162 }
163
164 Float_t AliPHOSPreprocessor::HG2LG(Int_t mod, Int_t X, Int_t Z, TFile* f)
165 {
166   //Calculates High gain to Low gain ratio 
167   //for crystal at the position (X,Z) in the PHOS module mod.
168   
169   char hname[128];
170   sprintf(hname,"%d_%d_%d",mod,X,Z);
171
172   TH1F* h1 = (TH1F*)f->Get(hname);
173   if(!h1) return 16.;
174
175   if(!h1->GetEntries()) return 16.;
176   if(h1->GetMaximum()<10.) return 16.;
177
178   Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
179   Double_t xmin = max - (h1->GetRMS()/3);
180   Double_t xmax = max + (h1->GetRMS()/2);
181   //       Double_t xmin = max - (h1->GetRMS());
182   //       Double_t xmax = max + (h1->GetRMS());
183
184   TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
185   gaus1->SetParNames("Constant","Mean","Sigma");
186   gaus1->SetParameter("Constant",h1->GetMaximum());
187   gaus1->SetParameter("Mean",max);
188   gaus1->SetParameter("Sigma",1.);
189   gaus1->SetLineColor(kBlue);
190   h1->Fit(gaus1,"LERQ+");
191
192   AliInfo(Form("%s: %d entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
193            h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),
194            gaus1->GetParameter("Mean"))); 
195
196   return gaus1->GetParameter("Mean");
197
198 }
199
200 Bool_t AliPHOSPreprocessor::FindBadChannelsEmc()
201 {
202   //Creates the bad channels map for PHOS EMC.
203
204   // The file fileName contains histograms which have been produced by DA2 detector algorithm.
205   // It is a responsibility of the SHUTTLE framework to form the fileName.
206
207   AliPHOSEmcBadChannelsMap badMap;
208
209   TList* list = GetFileSources(kDAQ, "BAD_CHANNELS");
210
211   if(!list) {
212     Log("Sources list for BAD_CHANNELS not found, exit.");
213     return kFALSE;
214   }
215
216   if(!list->GetEntries()) {
217     Log(Form("Got empty sources list. It seems DA2 did not produce any files!"));
218     return kFALSE;
219   }
220
221   TIter iter(list);
222   TObjString *source;
223   char hnam[80];
224   TH1F* h1=0;
225
226   const Float_t fQualityCut = 1.;
227   Int_t nGoods[5] = {0,0,0,0,0};
228   
229   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
230
231     AliInfo(Form("found source %s", source->String().Data()));
232
233     TString fileName = GetFile(kDAQ, "BAD_CHANNELS", source->GetName());
234     AliInfo(Form("Got filename: %s",fileName.Data()));
235
236     TFile f(fileName);
237
238     if(!f.IsOpen()) {
239       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
240       return kFALSE;
241     }
242     
243     Log(Form("Begin check for bad channels."));
244     
245     for(Int_t mod=0; mod<5; mod++) {
246       for(Int_t iX=0; iX<64; iX++) {
247         for(Int_t iZ=0; iZ<56; iZ++) {
248
249           sprintf(hnam,"%d_%d_%d_%d",mod,iX,iZ,1); // high gain 
250           h1 = (TH1F*)f.Get(hnam);
251
252           if(h1) {
253             Double_t mean = h1->GetMean();
254             
255             if(mean)
256               Log(Form("iX=%d iZ=%d gain=%d   mean=%.3f\n",iX,iZ,1,mean));
257
258             if( mean>0 && mean<fQualityCut ) { 
259               nGoods[mod]++; 
260             }
261             else
262               badMap.SetBadChannel(mod+1,iZ+1,iX+1); //module, col,row
263           }
264
265         }
266       }
267
268       if(nGoods[mod])
269         Log(Form("Module %d: %d good channels.",mod,nGoods[mod]));
270     }
271
272     
273   } // end of loop over sources
274   
275   // Store the bad channels map.
276   
277   AliCDBMetaData md;
278   md.SetResponsible("Boris Polishchuk");
279
280   // Data valid from current run fRun until being updated (validityInfinite=kTRUE)
281   Bool_t result = Store("Calib", "EmcBadChannels", &badMap, &md, fRun, kTRUE);
282   return result;
283
284 }
285
286 Bool_t AliPHOSPreprocessor::CalibrateEmc()
287 {
288   // I.  Calculates the set of calibration coefficients to equalyze the mean energies deposited at high gain.
289   // II. Extracts High_Gain/Low_Gain ratio for each channel.
290
291   // The file fileName contains histograms which have been produced by DA1 detector algorithm.
292   // It is a responsibility of the SHUTTLE framework to form the fileName.
293
294   gRandom->SetSeed(0); //the seed is set to the current  machine clock!
295   AliPHOSEmcCalibData calibData;
296   
297   TList* list = GetFileSources(kDAQ, "AMPLITUDES");
298   
299   if(!list) {
300     Log("Sources list not found, exit.");
301     return 1;
302   }
303
304   if(!list->GetEntries()) {
305     Log(Form("Got empty sources list. It seems DA1 did not produce any files!"));
306     return kFALSE;
307   }
308
309   TIter iter(list);
310   TObjString *source;
311   
312   while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
313     AliInfo(Form("found source %s", source->String().Data()));
314
315     TString fileName = GetFile(kDAQ, "AMPLITUDES", source->GetName());
316     AliInfo(Form("Got filename: %s",fileName.Data()));
317
318     TFile f(fileName);
319
320     if(!f.IsOpen()) {
321       Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
322       return 1;
323     }
324     
325     const Int_t nMod=5; // 1:5 modules
326     const Int_t nCol=56; //1:56 columns in each module
327     const Int_t nRow=64; //1:64 rows in each module
328
329     Double_t coeff;
330     char hnam[80];
331     TH2F* h2=0;
332     TH1D* h1=0;
333     
334     //Get the reference histogram
335     //(author: Gustavo Conesa Balbastre)
336
337     TList * keylist = f.GetListOfKeys();
338     Int_t nkeys   = f.GetNkeys();
339     Bool_t ok = kFALSE;
340     TKey  *key;
341     Int_t ikey = 0;
342     Int_t counter = 0;
343     TH1D* hRef = 0;
344
345     //Check if the file contains any histogram
346     
347     if(nkeys< 2){
348       Log(Form("Not enough histograms (%d) for calibration.",nkeys));
349       return 1;
350     }
351     
352     while(!ok){
353       ikey = gRandom->Integer(nkeys);
354       key = (TKey*)keylist->At(ikey);
355       TObject* obj = f.Get(key->GetName());
356       TString cname(obj->ClassName());
357       if(cname == "TH2F") {
358         h2 = (TH2F*)obj;
359         TString htitl = h2->GetTitle();
360         if(htitl.Contains("and gain 1")) {
361           hRef = h2->ProjectionX();
362           hRef->GetXaxis()->SetRange(10,1000); // to cut off saturation peak and noise
363           // Check if the reference histogram has too little statistics
364           if(hRef->GetMean() && hRef->GetEntries()>2) ok=kTRUE;
365         }
366       }
367       
368       counter++;
369       
370       if(!ok && counter > nkeys){
371         Log("No histogram with enough statistics for reference. Exit.");
372         return 1;
373       }
374     }
375     
376     Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
377              hRef->GetName(),hRef->GetEntries(),
378              hRef->GetMean(),hRef->GetRMS()));
379
380     Double_t refMean=hRef->GetMean();
381     
382     // Calculates relative calibration coefficients for all non-zero channels
383     
384     for(Int_t mod=0; mod<nMod; mod++) {
385       for(Int_t col=0; col<nCol; col++) {
386         for(Int_t row=0; row<nRow; row++) {
387
388           //High Gain to Low Gain ratio
389           Float_t ratio = HG2LG(mod,row,col,&f);
390           calibData.SetHighLowRatioEmc(mod+1,col+1,row+1,ratio);
391           
392           sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
393           h2 = (TH2F*)f.Get(hnam);
394
395           //TODO: dead channels exclusion!
396           if(h2) {
397             h1 = h2->ProjectionX();
398             h1->GetXaxis()->SetRange(10,1000); //to cut off saturation peak and noise
399             coeff = h1->GetMean()/refMean;
400             if(coeff>0)
401               calibData.SetADCchannelEmc(mod+1,col+1,row+1,1./coeff);
402             else 
403               calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.);
404             AliInfo(Form("mod %d col %d row %d  coeff %f\n",mod,col,row,coeff));
405           }
406           else
407             calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.); 
408         }
409       }
410     }
411     
412     f.Close();
413   }
414   
415   //Store EMC calibration data
416   
417   AliCDBMetaData emcMetaData;
418   Bool_t result = Store("Calib", "EmcGainPedestals", &calibData, &emcMetaData, 0, kTRUE); // valid from 0 to infinity);
419   return result;
420
421 }
422
423      
424