1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
23 // Author: Boris Polichtchouk, 4 October 2006
24 ///////////////////////////////////////////////////////////////////////////////
26 #include "AliPHOSPreprocessor.h"
28 #include "AliCDBMetaData.h"
29 #include "AliPHOSEmcCalibData.h"
38 #include "TObjString.h"
42 #include "AliPHOSEmcBadChannelsMap.h"
44 ClassImp(AliPHOSPreprocessor)
46 Double_t rgaus(Double_t *x, Double_t *par)
48 Double_t gaus = par[0] * TMath::Exp( -(x[0]-par[1])*(x[0]-par[1]) /
53 //_______________________________________________________________________________________
54 AliPHOSPreprocessor::AliPHOSPreprocessor() :
55 AliPreprocessor("PHS",0)
60 //_______________________________________________________________________________________
61 AliPHOSPreprocessor::AliPHOSPreprocessor(AliShuttleInterface* shuttle):
62 AliPreprocessor("PHS",shuttle)
67 //_______________________________________________________________________________________
68 UInt_t AliPHOSPreprocessor::Process(TMap* /*valueSet*/)
70 // process data retrieved by the Shuttle
72 // The fileName with the histograms which have been produced by
73 // AliPHOSCalibHistoProducer.
74 // It is a responsibility of the SHUTTLE framework to form the fileName
76 TString runType = GetRunType();
77 Log(Form("Run type: %s",runType.Data()));
80 Bool_t ledOK = ProcessLEDRun();
86 gRandom->SetSeed(0); //the seed is set to the current machine clock!
87 AliPHOSEmcCalibData calibData;
89 TList* list = GetFileSources(kDAQ, "AMPLITUDES");
91 Log("Sources list not found, exit.");
98 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
99 AliInfo(Form("found source %s", source->String().Data()));
101 TString fileName = GetFile(kDAQ, "AMPLITUDES", source->GetName());
102 AliInfo(Form("Got filename: %s",fileName.Data()));
107 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
111 const Int_t nMod=5; // 1:5 modules
112 const Int_t nCol=56; //1:56 columns in each module
113 const Int_t nRow=64; //1:64 rows in each module
120 //Get the reference histogram
121 //(author: Gustavo Conesa Balbastre)
123 TList * keylist = f.GetListOfKeys();
124 Int_t nkeys = f.GetNkeys();
131 //Check if the file contains any histogram
134 Log(Form("Not enough histograms (%d) for calibration.",nkeys));
139 ikey = gRandom->Integer(nkeys);
140 key = (TKey*)keylist->At(ikey);
141 TObject* obj = f.Get(key->GetName());
142 TString cname(obj->ClassName());
143 if(cname == "TH2F") {
145 TString htitl = h2->GetTitle();
146 if(htitl.Contains("and gain 1")) {
147 hRef = h2->ProjectionX();
148 hRef->SetBins(1000,0.,1000.); // to cut off saturation peak at 1023
149 // Check if the reference histogram has too little statistics
150 if(hRef->GetEntries()>2) ok=kTRUE;
151 if(!ok && counter > nkeys){
152 Log("No histogram with enough statistics for reference.");
160 Log(Form("reference histogram %s, %.1f entries, mean=%.3f, rms=%.3f.",
161 hRef->GetName(),hRef->GetEntries(),
162 hRef->GetMean(),hRef->GetRMS()));
164 Double_t refMean=hRef->GetMean();
166 // Calculates relative calibration coefficients for all non-zero channels
168 for(Int_t mod=0; mod<nMod; mod++) {
169 for(Int_t col=0; col<nCol; col++) {
170 for(Int_t row=0; row<nRow; row++) {
172 //High Gain to Low Gain ratio
173 Float_t ratio = HG2LG(mod,row,col,&f);
174 calibData.SetHighLowRatioEmc(mod+1,col+1,row+1,ratio);
176 sprintf(hnam,"%d_%d_%d_1",mod,row,col); // high gain!
177 h2 = (TH2F*)f.Get(hnam);
179 //TODO: dead channels exclusion!
181 h1 = h2->ProjectionX();
182 h1->SetBins(1000,0.,1000.); // to cut off saturation peak at 1023
183 coeff = h1->GetMean()/refMean;
185 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1./coeff);
187 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.);
188 AliInfo(Form("mod %d col %d row %d coeff %f\n",mod,col,row,coeff));
191 calibData.SetADCchannelEmc(mod+1,col+1,row+1,1.);
199 //Store EMC calibration data
201 AliCDBMetaData emcMetaData;
202 Bool_t emcOK = Store("Calib", "EmcGainPedestals", &calibData, &emcMetaData);
211 Bool_t AliPHOSPreprocessor::ProcessLEDRun()
213 //Process LED run, fill bad channels map.
215 AliPHOSEmcBadChannelsMap badMap;
217 TList* list = GetFileSources(kDAQ, "LED");
219 Log("Sources list for LED run not found, exit.");
228 while ((source = dynamic_cast<TObjString *> (iter.Next()))) {
230 AliInfo(Form("found source %s", source->String().Data()));
232 TString fileName = GetFile(kDAQ, "LED", source->GetName());
233 AliInfo(Form("Got filename: %s",fileName.Data()));
238 Log(Form("File %s is not opened, something goes wrong!",fileName.Data()));
242 const Int_t nMod=5; // 1:5 modules
243 const Int_t nCol=56; //1:56 columns in each module
244 const Int_t nRow=64; //1:64 rows in each module
246 // Check for dead channels
247 Log(Form("Begin check for dead channels."));
249 for(Int_t mod=0; mod<nMod; mod++) {
250 for(Int_t col=0; col<nCol; col++) {
251 for(Int_t row=0; row<nRow; row++) {
252 sprintf(hnam,"mod%dcol%drow%d",mod,col,row);
253 histo = (TH1F*)f.Get(hnam);
255 if (histo->GetMean()<1) {
256 Log(Form("Channel: [%d,%d,%d] seems dead, <E>=%.1f.",mod,col,row,histo->GetMean()));
257 badMap.SetBadChannel(mod,col,row);
265 //Store bad channels map
266 AliCDBMetaData badMapMetaData;
268 //Bad channels data valid from current run fRun until updated (validityInfinite=kTRUE)
269 Bool_t result = Store("Calib", "EmcBadChannels", &badMap, &badMapMetaData, fRun, kTRUE);
274 Float_t AliPHOSPreprocessor::HG2LG(Int_t mod, Int_t X, Int_t Z, TFile* f)
276 //Calculates High gain to Low gain ratio
277 //for crystal at the position (X,Z) in the PHOS module mod.
280 sprintf(hname,"%d_%d_%d",mod,X,Z);
282 TH1F* h1 = (TH1F*)f->Get(hname);
285 if(!h1->GetEntries()) return 16.;
286 if(h1->GetMaximum()<10.) return 16.;
288 Double_t max = h1->GetBinCenter(h1->GetMaximumBin()); // peak
289 Double_t xmin = max - (h1->GetRMS()/3);
290 Double_t xmax = max + (h1->GetRMS()/2);
291 // Double_t xmin = max - (h1->GetRMS());
292 // Double_t xmax = max + (h1->GetRMS());
294 TF1* gaus1 = new TF1("gaus1",rgaus,xmin,xmax,3);
295 gaus1->SetParNames("Constant","Mean","Sigma");
296 gaus1->SetParameter("Constant",h1->GetMaximum());
297 gaus1->SetParameter("Mean",max);
298 gaus1->SetParameter("Sigma",1.);
299 gaus1->SetLineColor(kBlue);
300 h1->Fit(gaus1,"LERQ+");
302 Log(Form("\t%s: %d entries, mean=%.3f, peak=%.3f, rms= %.3f. HG/LG = %.3f\n",
303 h1->GetTitle(),h1->GetEntries(),h1->GetMean(),max,h1->GetRMS(),
304 gaus1->GetParameter("Mean")));
306 return gaus1->GetParameter("Mean");