]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDPedestalDA.cxx
Fixing minor typo.
[u/mrichter/AliRoot.git] / FMD / AliFMDPedestalDA.cxx
CommitLineData
3bd993ba 1/**************************************************************************
2 * Copyright(c) 2004, 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/** @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.
20*/
21//
e9c06036 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.
3bd993ba 28//
29
30#include "AliFMDPedestalDA.h"
31#include "iostream"
32#include "fstream"
33#include "AliLog.h"
34#include "TF1.h"
35
36//_____________________________________________________________________
37ClassImp(AliFMDPedestalDA)
38
39//_____________________________________________________________________
40AliFMDPedestalDA::AliFMDPedestalDA() : AliFMDBaseDA()
41{
42 fOutputFile.open("peds.csv");
43
44}
45
46//_____________________________________________________________________
47AliFMDPedestalDA::AliFMDPedestalDA(const AliFMDPedestalDA & pedDA) :
48 AliFMDBaseDA(pedDA)
49{
50
51}
52
53//_____________________________________________________________________
e9c06036 54AliFMDPedestalDA::~AliFMDPedestalDA()
55{
3bd993ba 56}
57
58//_____________________________________________________________________
e9c06036 59void AliFMDPedestalDA::Init()
60{
3bd993ba 61 SetRequiredEvents(1000);
62}
63
64//_____________________________________________________________________
65void AliFMDPedestalDA::AddChannelContainer(TObjArray* sectorArray,
66 UShort_t det,
e9c06036 67 Char_t ring,
3bd993ba 68 UShort_t sec,
e9c06036 69 UShort_t strip)
70{
71#ifdef USE_SAMPLES
72 AliFMDParameters* pars = AliFMDParameters::Instance();
73 Int_t samples = pars->GetSampleRate(det, ring, sec, strip);
74 TObjArray* sampleArray = new TObjArray(samples);
75 for (size_t sample = 0; sample < samples; sample++) {
76 TH1S* hSample = new TH1S(Form("FMD%d%c[%02d,03%d]_%d",
77 det,ring,sec,strip,sample),
78 Form("FMD%d%c[%02d,%03%d]_%d",
79 det,ring,sec,strip),
80 1024,-.5,1023.5);
81 hSample->SetXTitle("ADC");
82 hSample->SetYTitle("Events");
83 sampleArray->AddAt(hSample, sample);
84 }
85 sectorArray->AddAtAndExpand(sampleArray, strip);
86#else
3bd993ba 87 TH1S* hChannel = new TH1S(Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
88 Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
e9c06036 89 1024,-.5,1023.5);
3bd993ba 90
91 hChannel->SetDirectory(0);
92 sectorArray->AddAtAndExpand(hChannel,strip);
e9c06036 93#endif
94
3bd993ba 95}
96
97//_____________________________________________________________________
e9c06036 98void AliFMDPedestalDA::FillChannels(AliFMDDigit* digit)
99{
3bd993ba 100 UShort_t det = digit->Detector();
101 Char_t ring = digit->Ring();
102 UShort_t sec = digit->Sector();
103 UShort_t strip = digit->Strip();
e9c06036 104
105#ifdef USE_SAMPLES
106 AliFMDParameters* pars = AliFMDParameters::Instance();
107 Int_t samples = pars->GetSampleRate(det, ring, sec, strip);
108 for (size_t sample = 0; sample < samples; sample++) {
109 TH1S* hSample = GetChannel(det, ring, sec, strip, sample);
110 hSample->Fill(digit->Count(sample));
111 }
112#else
113 TH1S* hChannel = GetChannel(det, ring, sec, strip);
3bd993ba 114 hChannel->Fill(digit->Counts());
e9c06036 115#endif
3bd993ba 116}
117
118//_____________________________________________________________________
119void AliFMDPedestalDA::Analyse(UShort_t det,
e9c06036 120 Char_t ring,
3bd993ba 121 UShort_t sec,
122 UShort_t strip) {
123
124 TH1S* hChannel = GetChannel(det, ring, sec, strip);
80fdb9f3 125 if(hChannel->GetEntries() == 0) {
e9c06036 126 // AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
127 // det,ring,sec,strip));
3bd993ba 128 return;
129 }
80fdb9f3 130
e9c06036 131 AliDebug(50, Form("Fitting FMD%d%c_%d_%d with %d entries",det,ring,sec,strip,
132 hChannel->GetEntries()));
3bd993ba 133
134 Float_t mean = hChannel->GetMean();
135 Float_t rms = hChannel->GetRMS();
136
80fdb9f3 137 hChannel->GetXaxis()->SetRangeUser(mean-4*rms,mean+4*rms);
3bd993ba 138
139 mean = hChannel->GetMean();
140 rms = hChannel->GetRMS();
141
e9c06036 142 hChannel->Fit("gaus","Q0+","",mean-5*rms,mean+5*rms);
3bd993ba 143 TF1* fitFunc = hChannel->GetFunction("gaus");
144
145 UInt_t ddl, board, chip, channel;
146
e9c06036 147 UShort_t relStrip = strip%128;
148 AliFMDParameters* pars = AliFMDParameters::Instance();
3bd993ba 149 pars->Detector2Hardware(det,ring,sec,strip,ddl,board,chip,channel);
e9c06036 150
3bd993ba 151 Float_t chi2ndf = 0;
152 if(fitFunc->GetNDF())
153 chi2ndf = fitFunc->GetChisquare() / fitFunc->GetNDF();
e9c06036 154
3bd993ba 155 ddl = ddl + kBaseDDL;
156
157 fOutputFile << ddl << ','
158 << board << ','
159 << chip << ','
160 << channel << ','
161 << relStrip << ','
162 << 0 << ','
163 << 0 << ','
164 << mean << ','
165 << rms << ','
166 << fitFunc->GetParameter(1) << ','
167 << fitFunc->GetParameter(2) << ','
168 << chi2ndf <<"\n";
169
170 if(fSaveHistograms) {
e9c06036 171 gDirectory->cd(GetSectorPath(det, ring, sec, kTRUE));
172 TH1F* sumPed = dynamic_cast<TH1F*>(gDirectory->Get("Pedestals"));
173 TH1F* sumNoise = dynamic_cast<TH1F*>(gDirectory->Get("Noise"));
174 Int_t nStr = (ring == 'I' ? 512 : 256);
175 if (!sumPed) {
176 sumPed = new TH1F("Pedestals",
177 Form("Summary of pedestals in FMD%d%c[%02d]",
178 det, ring, sec),
179 nStr, -.5, nStr-.5);
180 sumPed->SetXTitle("Strip");
181 sumPed->SetYTitle("Pedestal [ADC]");
182 sumPed->SetDirectory(gDirectory);
183 }
184 if (!sumNoise) {
185 sumNoise = new TH1F("Noise",
186 Form("Summary of noise in FMD%d%c[%02d]",
187 det, ring, sec),
188 nStr, -.5, nStr-.5);
189 sumNoise->SetXTitle("Strip");
190 sumNoise->SetYTitle("Noise [ADC]");
191 sumNoise->SetDirectory(gDirectory);
192 }
193 sumPed->SetBinContent(strip+1, mean);
194 sumPed->SetBinError(strip+1, rms);
195 sumNoise->SetBinContent(strip+1, rms);
196
197 gDirectory->cd(GetStripPath(det, ring, sec, strip, kTRUE));
198 hChannel->GetXaxis()->SetRange(-.5,1023.5);
199 // hChannel->SetDirectory(gDirectory);
3bd993ba 200 hChannel->Write();
e9c06036 201 }
3bd993ba 202}
203
204//_____________________________________________________________________
e9c06036 205void AliFMDPedestalDA::WriteHeaderToFile()
206{
80fdb9f3 207 AliFMDParameters* pars = AliFMDParameters::Instance();
208 fOutputFile.write(Form("# %s \n",pars->GetPedestalShuttleID()),13);
e9c06036 209 fOutputFile.write("# Rcu, "
210 "Board, "
211 "Chip, "
212 "Channel, "
213 "Strip, "
214 "Sample, "
215 "TimeBin, "
216 "Pedestal, "
217 "Noise, "
218 "Mu, "
219 "Sigma, "
220 "Chi2/NDF \n",91);
3bd993ba 221}
222
223//_____________________________________________________________________
e9c06036 224TH1S* AliFMDPedestalDA::GetChannel(UShort_t det,
225 Char_t ring,
226 UShort_t sec,
227 UShort_t strip)
228{
229 UShort_t iring = (ring == 'O' ? 0 : 1);
3bd993ba 230 TObjArray* detArray = static_cast<TObjArray*>(fDetectorArray.At(det));
e9c06036 231 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(iring));
3bd993ba 232 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
e9c06036 233#ifdef USE_SAMPLES
234 AliFMDParameters* pars = AliFMDParameters::Instance();
235 Int_t samples = pars->GetSampleRate(det, ring, sec, strip);
236 TObjArray* sampleArray = static_cast<TObjArray*>(secArray->At(strip));
237 TH1S* hSample = static_cast<TH1S*>(sampleArray->At(sample));
238 return hSample;
239#else
240 TH1S* hChannel = static_cast<TH1S*>(secArray->At(strip));
3bd993ba 241 return hChannel;
e9c06036 242#endif
3bd993ba 243}
244//_____________________________________________________________________
245//
246//EOF
247//