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"
38 //_____________________________________________________________________
39 ClassImp(AliFMDPedestalDA)
41 //_____________________________________________________________________
42 AliFMDPedestalDA::AliFMDPedestalDA() : AliFMDBaseDA(),
44 fPedSummary("PedestalSummary","pedestals",51200,0,51200),
45 fNoiseSummary("NoiseSummary","noise",51200,0,51200),
50 fOutputFile.open("peds.csv");
51 fZSfileFMD1.open("ddl3072.csv");
52 fZSfileFMD2.open("ddl3073.csv");
53 fZSfileFMD3.open("ddl3074.csv");
56 //_____________________________________________________________________
57 AliFMDPedestalDA::AliFMDPedestalDA(const AliFMDPedestalDA & pedDA) :
60 fPedSummary("PedestalSummary","pedestals",51200,0,51200),
61 fNoiseSummary("NoiseSummary","noise",51200,0,51200),
69 //_____________________________________________________________________
70 AliFMDPedestalDA::~AliFMDPedestalDA()
74 //_____________________________________________________________________
75 void AliFMDPedestalDA::Init()
77 SetRequiredEvents(1000);
82 //_____________________________________________________________________
83 void AliFMDPedestalDA::AddChannelContainer(TObjArray* sectorArray,
89 AliFMDParameters* pars = AliFMDParameters::Instance();
90 UInt_t samples = pars->GetSampleRate(det, ring, sec, strip);
91 TObjArray* sampleArray = new TObjArray(samples);
92 sampleArray->SetOwner();
93 for (UInt_t sample = 0; sample < samples; sample++) {
94 TH1S* hSample = new TH1S(Form("FMD%d%c[%02d,03%d]_%d",
95 det,ring,sec,strip,sample),
96 Form("FMD%d%c[%02d,%03%d]_%d",
99 hSample->SetXTitle("ADC");
100 hSample->SetYTitle("Events");
101 hSample->SetDirectory(0);
102 sampleArray->AddAt(hSample, sample);
104 sectorArray->AddAtAndExpand(sampleArray, strip);
109 //_____________________________________________________________________
110 void AliFMDPedestalDA::FillChannels(AliFMDDigit* digit)
112 UShort_t det = digit->Detector();
113 Char_t ring = digit->Ring();
114 UShort_t sec = digit->Sector();
115 UShort_t strip = digit->Strip();
117 AliFMDParameters* pars = AliFMDParameters::Instance();
118 UInt_t samples = pars->GetSampleRate(det, ring, sec, strip);
119 for (UInt_t sample = 0; sample < samples; sample++) {
120 TH1S* hSample = GetChannel(det, ring, sec, strip, sample);
121 hSample->Fill(digit->Count(sample));
126 //_____________________________________________________________________
127 void AliFMDPedestalDA::Analyse(UShort_t det,
132 AliFMDParameters* pars = AliFMDParameters::Instance();
133 UInt_t samples = pars->GetSampleRate(det, ring, sec, strip);
134 for (UShort_t sample = 0; sample < samples; sample++) {
136 TH1S* hChannel = GetChannel(det, ring, sec, strip,sample);
137 if(hChannel->GetEntries() == 0) {
138 //AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
139 // det,ring,sec,strip));
143 AliDebug(50, Form("Fitting FMD%d%c_%d_%d with %d entries",det,ring,sec,strip,
144 hChannel->GetEntries()));
145 TF1 fitFunc("fitFunc","gausn",0,300);
146 fitFunc.SetParameters(100,100,1);
147 hChannel->Fit("fitFunc","Q0+","",10,200);
149 Float_t mean = hChannel->GetMean();
150 Float_t rms = hChannel->GetRMS();
154 hChannel->GetXaxis()->SetRangeUser(mean-5*rms,mean+5*rms);
156 mean = hChannel->GetMean();
157 rms = hChannel->GetRMS();
160 UShort_t ddl, board, altro, channel;
163 pars->Detector2Hardware(det,ring,sec,strip,sample,ddl,board,altro,channel,timebin);
167 fZSfileFMD1 << board << ',' << altro << ',' << channel << ',' << timebin << ','
168 << mean << ',' << rms << "\n"; break;
170 fZSfileFMD2 << board << ',' << altro << ',' << channel << ',' << timebin << ','
171 << mean << ',' << rms << "\n"; break;
173 fZSfileFMD3 << board << ',' << altro << ',' << channel << ',' << timebin << ','
174 << mean << ',' << rms << "\n"; break;
176 AliWarning("Unknown sample!"); break;
184 chi2ndf = fitFunc.GetChisquare() / fitFunc.GetNDF();
186 Int_t sampleToWrite = 2;
188 if(pars->GetSampleRate(det,ring,sec,strip)==2)
191 if(pars->GetSampleRate(det,ring,sec,strip)<2)
194 if(sample==sampleToWrite) {
196 fOutputFile << det << ','
202 << fitFunc.GetParameter(1) << ','
203 << fitFunc.GetParameter(2) << ','
206 if(fSaveHistograms ) {
207 gDirectory->cd(GetSectorPath(det, ring, sec, kTRUE));
208 TH1F* sumPed = dynamic_cast<TH1F*>(gDirectory->Get("Pedestals"));
209 TH1F* sumNoise = dynamic_cast<TH1F*>(gDirectory->Get("Noise"));
210 Int_t nStr = (ring == 'I' ? 512 : 256);
212 sumPed = new TH1F("Pedestals",
213 Form("Summary of pedestals in FMD%d%c[%02d]",
216 sumPed->SetXTitle("Strip");
217 sumPed->SetYTitle("Pedestal [ADC]");
218 sumPed->SetDirectory(gDirectory);
221 sumNoise = new TH1F("Noise",
222 Form("Summary of noise in FMD%d%c[%02d]",
225 sumNoise->SetXTitle("Strip");
226 sumNoise->SetYTitle("Noise [ADC]");
228 sumNoise->SetDirectory(gDirectory);
230 sumPed->SetBinContent(strip+1, mean);
231 sumPed->SetBinError(strip+1, rms);
232 sumNoise->SetBinContent(strip+1, rms);
234 if(sumNoise->GetEntries() == nStr)
235 sumNoise->Write(sumNoise->GetName(),TObject::kOverwrite);
236 if(sumPed->GetEntries() == nStr)
237 sumPed->Write(sumPed->GetName(),TObject::kOverwrite);
239 fPedSummary.SetBinContent(fCurrentChannel,mean);
241 fNoiseSummary.SetBinContent(fCurrentChannel,rms);
244 gDirectory->cd(GetStripPath(det, ring, sec, strip, kTRUE));
245 hChannel->GetXaxis()->SetRange(1,1024);
253 //_____________________________________________________________________
254 void AliFMDPedestalDA::Terminate(TFile* diagFile)
256 if(fSaveHistograms) {
260 fNoiseSummary.Write();
263 if(fZSfileFMD1.is_open()) {
264 fZSfileFMD1.write("# EOF\n",6);
265 fZSfileFMD1.close(); }
266 if(fZSfileFMD2.is_open()) {
267 fZSfileFMD2.write("# EOF\n",6);
268 fZSfileFMD2.close(); }
269 if(fZSfileFMD3.is_open()) {
270 fZSfileFMD3.write("# EOF\n",6);
271 fZSfileFMD3.close(); }
275 //_____________________________________________________________________
276 void AliFMDPedestalDA::WriteHeaderToFile()
278 AliFMDParameters* pars = AliFMDParameters::Instance();
279 fOutputFile.write(Form("# %s \n",pars->GetPedestalShuttleID()),13);
280 fOutputFile.write("# Detector, "
289 fZSfileFMD1.write("# FMD 1 pedestals \n",19);
290 fZSfileFMD1.write("# board, "
296 fZSfileFMD2.write("# FMD 2 pedestals \n",19);
297 fZSfileFMD2.write("# board, "
303 fZSfileFMD3.write("# FMD 3 pedestals \n",19);
304 fZSfileFMD3.write("# board, "
313 //_____________________________________________________________________
314 TH1S* AliFMDPedestalDA::GetChannel(UShort_t det,
320 UShort_t iring = (ring == 'O' ? 0 : 1);
321 TObjArray* detArray = static_cast<TObjArray*>(fDetectorArray.At(det));
322 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(iring));
323 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
324 TObjArray* sampleArray = static_cast<TObjArray*>(secArray->At(strip));
325 TH1S* hSample = static_cast<TH1S*>(sampleArray->At(sample));
330 //_____________________________________________________________________