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"
37 //_____________________________________________________________________
38 ClassImp(AliFMDPedestalDA)
40 //_____________________________________________________________________
41 AliFMDPedestalDA::AliFMDPedestalDA() : AliFMDBaseDA(),
43 fPedSummary("PedestalSummary","pedestals",51200,0,51200),
44 fNoiseSummary("NoiseSummary","noise",51200,0,51200)
46 fOutputFile.open("peds.csv");
50 //_____________________________________________________________________
51 AliFMDPedestalDA::AliFMDPedestalDA(const AliFMDPedestalDA & pedDA) :
54 fPedSummary("PedestalSummary","pedestals",51200,0,51200),
55 fNoiseSummary("NoiseSummary","noise",51200,0,51200)
60 //_____________________________________________________________________
61 AliFMDPedestalDA::~AliFMDPedestalDA()
65 //_____________________________________________________________________
66 void AliFMDPedestalDA::Init()
68 SetRequiredEvents(1000);
71 //_____________________________________________________________________
72 void AliFMDPedestalDA::AddChannelContainer(TObjArray* sectorArray,
79 AliFMDParameters* pars = AliFMDParameters::Instance();
80 Int_t samples = pars->GetSampleRate(det, ring, sec, strip);
81 TObjArray* sampleArray = new TObjArray(samples);
82 for (size_t sample = 0; sample < samples; sample++) {
83 TH1S* hSample = new TH1S(Form("FMD%d%c[%02d,03%d]_%d",
84 det,ring,sec,strip,sample),
85 Form("FMD%d%c[%02d,%03%d]_%d",
88 hSample->SetXTitle("ADC");
89 hSample->SetYTitle("Events");
90 sampleArray->AddAt(hSample, sample);
92 sectorArray->AddAtAndExpand(sampleArray, strip);
94 TH1S* hChannel = new TH1S(Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
95 Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
98 hChannel->SetDirectory(0);
99 sectorArray->AddAtAndExpand(hChannel,strip);
104 //_____________________________________________________________________
105 void AliFMDPedestalDA::FillChannels(AliFMDDigit* digit)
107 UShort_t det = digit->Detector();
108 Char_t ring = digit->Ring();
109 UShort_t sec = digit->Sector();
110 UShort_t strip = digit->Strip();
113 AliFMDParameters* pars = AliFMDParameters::Instance();
114 Int_t samples = pars->GetSampleRate(det, ring, sec, strip);
115 for (size_t sample = 0; sample < samples; sample++) {
116 TH1S* hSample = GetChannel(det, ring, sec, strip, sample);
117 hSample->Fill(digit->Count(sample));
120 TH1S* hChannel = GetChannel(det, ring, sec, strip);
121 hChannel->Fill(digit->Counts());
125 //_____________________________________________________________________
126 void AliFMDPedestalDA::Analyse(UShort_t det,
131 TH1S* hChannel = GetChannel(det, ring, sec, strip);
132 if(hChannel->GetEntries() == 0) {
133 // AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
134 // det,ring,sec,strip));
138 AliDebug(50, Form("Fitting FMD%d%c_%d_%d with %d entries",det,ring,sec,strip,
139 hChannel->GetEntries()));
141 Float_t mean = hChannel->GetMean();
142 Float_t rms = hChannel->GetRMS();
144 hChannel->GetXaxis()->SetRangeUser(mean-4*rms,mean+4*rms);
146 mean = hChannel->GetMean();
147 rms = hChannel->GetRMS();
149 hChannel->Fit("gaus","Q0+","",mean-5*rms,mean+5*rms);
150 TF1* fitFunc = hChannel->GetFunction("gaus");
153 if(fitFunc->GetNDF())
154 chi2ndf = fitFunc->GetChisquare() / fitFunc->GetNDF();
156 fOutputFile << det << ','
162 << fitFunc->GetParameter(1) << ','
163 << fitFunc->GetParameter(2) << ','
166 if(fSaveHistograms) {
167 gDirectory->cd(GetSectorPath(det, ring, sec, kTRUE));
168 TH1F* sumPed = dynamic_cast<TH1F*>(gDirectory->Get("Pedestals"));
169 TH1F* sumNoise = dynamic_cast<TH1F*>(gDirectory->Get("Noise"));
170 Int_t nStr = (ring == 'I' ? 512 : 256);
172 sumPed = new TH1F("Pedestals",
173 Form("Summary of pedestals in FMD%d%c[%02d]",
176 sumPed->SetXTitle("Strip");
177 sumPed->SetYTitle("Pedestal [ADC]");
178 sumPed->SetDirectory(gDirectory);
181 sumNoise = new TH1F("Noise",
182 Form("Summary of noise in FMD%d%c[%02d]",
185 sumNoise->SetXTitle("Strip");
186 sumNoise->SetYTitle("Noise [ADC]");
188 sumNoise->SetDirectory(gDirectory);
190 sumPed->SetBinContent(strip+1, mean);
191 sumPed->SetBinError(strip+1, rms);
192 sumNoise->SetBinContent(strip+1, rms);
194 if(sumNoise->GetEntries() == nStr)
195 sumNoise->Write(sumNoise->GetName(),TObject::kOverwrite);
196 if(sumPed->GetEntries() == nStr)
197 sumPed->Write(sumPed->GetName(),TObject::kOverwrite);
199 fPedSummary.SetBinContent(fCurrentChannel,mean);
200 fNoiseSummary.SetBinContent(fCurrentChannel,rms);
203 gDirectory->cd(GetStripPath(det, ring, sec, strip, kTRUE));
204 hChannel->GetXaxis()->SetRange(1,1024);
210 //_____________________________________________________________________
211 void AliFMDPedestalDA::Terminate(TFile* diagFile)
216 fNoiseSummary.Write();
220 //_____________________________________________________________________
221 void AliFMDPedestalDA::WriteHeaderToFile()
223 AliFMDParameters* pars = AliFMDParameters::Instance();
224 fOutputFile.write(Form("# %s \n",pars->GetPedestalShuttleID()),13);
225 fOutputFile.write("# Detector, "
236 //_____________________________________________________________________
237 TH1S* AliFMDPedestalDA::GetChannel(UShort_t det,
242 UShort_t iring = (ring == 'O' ? 0 : 1);
243 TObjArray* detArray = static_cast<TObjArray*>(fDetectorArray.At(det));
244 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(iring));
245 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
247 AliFMDParameters* pars = AliFMDParameters::Instance();
248 Int_t samples = pars->GetSampleRate(det, ring, sec, strip);
249 TObjArray* sampleArray = static_cast<TObjArray*>(secArray->At(strip));
250 TH1S* hSample = static_cast<TH1S*>(sampleArray->At(sample));
253 TH1S* hChannel = static_cast<TH1S*>(secArray->At(strip));
257 //_____________________________________________________________________