1 /**************************************************************************
2 * Copyright(c) 2004, 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 **************************************************************************/
16 /** @file AliFMDPedestalDA.cxx
17 @author Hans Hjersing Dalsgaard <canute@nbi.dk>
18 @date Mon Mar 10 09:46:05 2008
19 @brief Derived class for the pedestal detector algorithm.
22 // This class implements the virtual functions of the AliFMDBaseDA
23 // class. The most important of these functions, FillChannels(..) and
24 // Analyse(...) collect and analyse the data of each channel. The
25 // resulting pedestal and noise values are written to a comma
26 // separated values (csv) file on the go. The csv files produced in
27 // this way are the basic input to the AliFMDPreprocessor.
30 #include "AliFMDPedestalDA.h"
31 #include "AliFMDAltroMapping.h"
40 //_____________________________________________________________________
41 ClassImp(AliFMDPedestalDA)
43 //_____________________________________________________________________
44 AliFMDPedestalDA::AliFMDPedestalDA() : AliFMDBaseDA(),
46 fPedSummary("PedestalSummary","pedestals",51200,0,51200),
47 fNoiseSummary("NoiseSummary","noise",51200,0,51200),
51 fMinTimebin(3 * 4 * 3 * 16), // 3 ddls, 4 FECs, 3 Altros, 16 channels
52 fMaxTimebin(3 * 4 * 3 * 16) // 3 ddls, 4 FECs, 3 Altros, 16 channels
54 // Default constructor
55 fOutputFile.open("peds.csv");
56 fZSfileFMD1.open("ddl3072.csv");
57 fZSfileFMD2.open("ddl3073.csv");
58 fZSfileFMD3.open("ddl3074.csv");
61 //_____________________________________________________________________
62 AliFMDPedestalDA::AliFMDPedestalDA(const AliFMDPedestalDA & pedDA) :
65 fPedSummary("PedestalSummary","pedestals",51200,0,51200),
66 fNoiseSummary("NoiseSummary","noise",51200,0,51200),
70 fMinTimebin(pedDA.fMinTimebin),
71 fMaxTimebin(pedDA.fMaxTimebin)
76 //_____________________________________________________________________
77 AliFMDPedestalDA::~AliFMDPedestalDA()
82 //_____________________________________________________________________
83 void AliFMDPedestalDA::Init()
86 SetRequiredEvents(1000);
87 fMinTimebin.Reset(1024);
88 fMaxTimebin.Reset(-1);
91 //_____________________________________________________________________
92 void AliFMDPedestalDA::AddChannelContainer(TObjArray* sectorArray,
98 // Add a channel to the containers.
101 // sectorArray Array of sectors
106 AliFMDParameters* pars = AliFMDParameters::Instance();
107 UInt_t samples = pars->GetSampleRate(det, ring, sec, strip);
108 TObjArray* sampleArray = new TObjArray(samples);
109 sampleArray->SetOwner();
110 for (UInt_t sample = 0; sample < samples; sample++) {
111 TH1S* hSample = new TH1S(Form("FMD%d%c[%02d,03%d]_%d",
112 det,ring,sec,strip,sample),
113 Form("FMD%d%c[%02d,%03%d]_%d",
116 hSample->SetXTitle("ADC");
117 hSample->SetYTitle("Events");
118 hSample->SetDirectory(0);
119 sampleArray->AddAt(hSample, sample);
121 sectorArray->AddAtAndExpand(sampleArray, strip);
124 //_____________________________________________________________________
125 void AliFMDPedestalDA::FillChannels(AliFMDDigit* digit)
127 // Fill ADC values from a digit into the corresponding histogram.
130 // digit Digit to fill ADC values for.
131 UShort_t det = digit->Detector();
132 Char_t ring = digit->Ring();
133 UShort_t sec = digit->Sector();
134 UShort_t strip = digit->Strip();
136 AliFMDParameters* pars = AliFMDParameters::Instance();
137 UInt_t samples = pars->GetSampleRate(det, ring, sec, strip);
138 for (UInt_t sample = 0; sample < samples; sample++) {
139 TH1S* hSample = GetChannel(det, ring, sec, strip, sample);
140 hSample->Fill(digit->Count(sample));
145 //_____________________________________________________________________
146 void AliFMDPedestalDA::Analyse(UShort_t det,
151 // Analyse a strip. That is, compute the mean and spread of the ADC
152 // spectra for all strips. Also output on files the values.
159 AliFMDParameters* pars = AliFMDParameters::Instance();
160 // Float_t factor = pars->GetPedestalFactor();
161 UInt_t samples = pars->GetSampleRate(det, ring, sec, strip);
162 for (UShort_t sample = 0; sample < samples; sample++) {
164 TH1S* hChannel = GetChannel(det, ring, sec, strip,sample);
165 if(hChannel->GetEntries() == 0) {
166 //AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
167 // det,ring,sec,strip));
171 AliDebug(50, Form("Fitting FMD%d%c_%d_%d with %d entries",
172 det,ring,sec,strip, hChannel->GetEntries()));
173 TF1 fitFunc("fitFunc","gausn",0,300);
174 fitFunc.SetParameters(100,100,1);
175 hChannel->Fit("fitFunc","Q0+","",10,200);
177 Float_t mean = hChannel->GetMean();
178 Float_t rms = hChannel->GetRMS();
182 hChannel->GetXaxis()->SetRangeUser(mean-5*rms,mean+5*rms);
184 mean = hChannel->GetMean();
185 rms = hChannel->GetRMS();
188 UShort_t ddl, board, altro, channel;
191 pars->Detector2Hardware(det,ring,sec,strip,sample,
192 ddl,board,altro,channel,timebin);
193 Int_t idx = HWIndex(ddl, board, altro, channel);
195 fMinTimebin[idx] = TMath::Min(Short_t(timebin), fMinTimebin[idx]);
196 fMaxTimebin[idx] = TMath::Max(Short_t(timebin+1), fMaxTimebin[idx]);
199 std::ostream* zsFile = 0;
201 case 1: zsFile = &fZSfileFMD1; break;
202 case 2: zsFile = &fZSfileFMD2; break;
203 case 3: zsFile = &fZSfileFMD3; break;
204 default: AliWarning("Unknown sample!"); break;
207 *zsFile << board << ','
218 chi2ndf = fitFunc.GetChisquare() / fitFunc.GetNDF();
221 Int_t sampleToWrite = 2;
222 if (samples == 2) sampleToWrite = 1;
223 else if (samples < 2) sampleToWrite = 0;
225 if(sample != sampleToWrite) continue;
228 fOutputFile << det << ','
234 << fitFunc.GetParameter(1) << ','
235 << fitFunc.GetParameter(2) << ','
238 if(fSaveHistograms ) {
239 gDirectory->cd(GetSectorPath(det, ring, sec, kTRUE));
240 TH1F* sumPed = dynamic_cast<TH1F*>(gDirectory->Get("Pedestals"));
241 TH1F* sumNoise = dynamic_cast<TH1F*>(gDirectory->Get("Noise"));
242 Int_t nStr = (ring == 'I' ? 512 : 256);
244 sumPed = new TH1F("Pedestals",
245 Form("Summary of pedestals in FMD%d%c[%02d]",
248 sumPed->SetXTitle("Strip");
249 sumPed->SetYTitle("Pedestal [ADC]");
250 sumPed->SetDirectory(gDirectory);
253 sumNoise = new TH1F("Noise",
254 Form("Summary of noise in FMD%d%c[%02d]",
257 sumNoise->SetXTitle("Strip");
258 sumNoise->SetYTitle("Noise [ADC]");
260 sumNoise->SetDirectory(gDirectory);
262 sumPed->SetBinContent(strip+1, mean);
263 sumPed->SetBinError(strip+1, rms);
264 sumNoise->SetBinContent(strip+1, rms);
266 if(sumNoise->GetEntries() == nStr)
267 sumNoise->Write(sumNoise->GetName(),TObject::kOverwrite);
268 if(sumPed->GetEntries() == nStr)
269 sumPed->Write(sumPed->GetName(),TObject::kOverwrite);
271 fPedSummary.SetBinContent(fCurrentChannel,mean);
273 fNoiseSummary.SetBinContent(fCurrentChannel,rms);
276 gDirectory->cd(GetStripPath(det, ring, sec, strip, kTRUE));
277 hChannel->GetXaxis()->SetRange(1,1024);
284 //_____________________________________________________________________
285 void AliFMDPedestalDA::Terminate(TFile* diagFile)
287 // Called at the end of a job. Fills in missing time-bins and
288 // closes output files
289 if(fSaveHistograms) {
293 fNoiseSummary.Write();
295 AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
296 for (Int_t i = 0; i < 3; i++) {
297 std::ofstream& out = (i == 0 ? fZSfileFMD1 :
298 i == 1 ? fZSfileFMD2 :
300 if (out.is_open() && fSeenDetectors[i]) {
301 FillinTimebins(out, map->Detector2DDL(i+1));
303 if (!fSeenDetectors[i]) {
304 TString n(Form("ddl%d.csv",3072+map->Detector2DDL(i+1)));
305 gSystem->Unlink(n.Data());
311 //_____________________________________________________________________
312 void AliFMDPedestalDA::FillinTimebins(std::ofstream& out, UShort_t /*ddl*/)
315 unsigned short boards[] = { 0x0, 0x1, 0x10, 0x11, 0xFFFF };
316 unsigned short* board = boards;
317 while ((*boards) != 0xFFFF) {
318 for (UShort_t altro = 0; altro < 3; altro++) {
319 for (UShort_t channel = 0; channel < 16; channel++) {
320 Int_t idx = HWIndex(ddl, *board, altro, channel);
322 AliWarning(Form("Invalid index for %4d/0x%02x/0x%x/0x%x: %d",
323 ddl, *board, altro, channel, idx));
326 Short_t min = fMinTimebin[idx];
327 Short_t max = fMaxTimebin[idx];
329 // Channel not seen at all.
330 if (min > 1023 || max < 0) continue;
332 out << "# Extra timebins for 0x" << std::hex
333 << board << ',' << altro << ',' << channel
334 << " got time-bins " << min << " to " << max-1
335 << std::dec << std::endl;
337 for (UShort_t t = 15; t < min; t++)
338 // Write a phony line
339 out << board << "," << altro << "," << channel << ","
340 << t << "," << 1023 << "," << 0 << std::endl;
342 for (UShort_t t = max; t < 1024; t++)
343 // Write a phony line
344 out << board << "," << altro << "," << channel << ","
345 << t << "," << 1023 << "," << 0 << std::endl;
349 // Write trailer, and close
351 out.write("# EOF\n", 6);
355 //_____________________________________________________________________
356 void AliFMDPedestalDA::WriteHeaderToFile()
358 // Write headers to output files
359 AliFMDParameters* pars = AliFMDParameters::Instance();
360 fOutputFile.write(Form("# %s \n",pars->GetPedestalShuttleID()),13);
361 fOutputFile.write("# Detector, "
371 std::ostream* zss[] = { &fZSfileFMD1, &fZSfileFMD2, &fZSfileFMD3, 0 };
372 for (size_t i = 0; i < 3; i++)
373 *(zss[i]) << "# FMD 1 pedestals \n"
382 //_____________________________________________________________________
383 TH1S* AliFMDPedestalDA::GetChannel(UShort_t det,
389 // Get the histogram corresponding to a strip sample.
399 // ADC spectra of a strip.
400 UShort_t iring = (ring == 'O' ? 0 : 1);
401 TObjArray* detArray = static_cast<TObjArray*>(fDetectorArray.At(det));
402 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(iring));
403 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
404 TObjArray* sampleArray = static_cast<TObjArray*>(secArray->At(strip));
405 TH1S* hSample = static_cast<TH1S*>(sampleArray->At(sample));
410 //_____________________________________________________________________