]> git.uio.no Git - u/mrichter/AliRoot.git/blame - FMD/AliFMDGainDA.cxx
Removed a cout statement
[u/mrichter/AliRoot.git] / FMD / AliFMDGainDA.cxx
CommitLineData
f631f26d 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 AliFMDGainDA.cxx
17 @author Hans Hjersing Dalsgaard <canute@nbi.dk>
18 @date Mon Mar 13 13:46:05 2008
19 @brief Derived class for the pulse gain detector algorithm.
20*/
cf291b42 21//____________________________________________________________________
22//
23// This class contains the implementation of the gain detector
24// algorithms (DA) for the FMD. The data is collected in histograms
25// that are reset for each pulse length after the mean and standard
26// deviation are put into a TGraphErrors object. After a certain
27// number of pulses (usually 8) the graph is fitted to a straight
28// line. The gain is then slope of this line as it combines the known
29// pulse and the response of the detector.
30//
f631f26d 31#include "AliFMDGainDA.h"
32#include "iostream"
33#include "fstream"
34#include "AliLog.h"
35#include "TF1.h"
36#include "TH1.h"
37#include "TMath.h"
38#include "TGraphErrors.h"
39#include "AliFMDParameters.h"
40
41//_____________________________________________________________________
42ClassImp(AliFMDGainDA)
cf291b42 43#if 0 // Do not delete - here to let Emacs indent properly
44;
45#endif
f631f26d 46
47//_____________________________________________________________________
cf291b42 48AliFMDGainDA::AliFMDGainDA()
49 : AliFMDBaseDA(),
50 fGainArray(),
427e8f99 51 //fPulseSize(32),
cf291b42 52 fHighPulse(256),
427e8f99 53 //fPulseLength(100),
54 fEventsPerChannel(16),
55 fCurrentPulse(16),
56 fCurrentChannel(16),
3bae5d02 57 fNumberOfStripsPerChip(128),
58 fSummaryGains("GainsSummary","Summary of gains",51200,0,51200),
f14ede67 59 fCurrentSummaryStrip(1)
f631f26d 60{
427e8f99 61 fCurrentPulse.Reset(0);
62 fCurrentChannel.Reset(0);
f631f26d 63 fOutputFile.open("gains.csv");
cf291b42 64 fGainArray.SetOwner();
f631f26d 65}
66
67//_____________________________________________________________________
cf291b42 68AliFMDGainDA::AliFMDGainDA(const AliFMDGainDA & gainDA)
69 : AliFMDBaseDA(gainDA),
70 fGainArray(gainDA.fGainArray),
427e8f99 71 //fPulseSize(gainDA.fPulseSize),
cf291b42 72 fHighPulse(gainDA.fHighPulse),
427e8f99 73 //fPulseLength(gainDA.fPulseLength),
cf291b42 74 fEventsPerChannel(gainDA.fEventsPerChannel),
75 fCurrentPulse(gainDA.fCurrentPulse),
76 fCurrentChannel(gainDA.fCurrentChannel),
3bae5d02 77 fNumberOfStripsPerChip(gainDA.fNumberOfStripsPerChip),
78 fSummaryGains(gainDA.fSummaryGains),
79 fCurrentSummaryStrip(gainDA.fCurrentSummaryStrip)
cf291b42 80{
427e8f99 81 fCurrentPulse.Reset(0);
82 fCurrentChannel.Reset(0);
f631f26d 83}
84
85//_____________________________________________________________________
cf291b42 86AliFMDGainDA::~AliFMDGainDA()
87{
f631f26d 88}
89
90//_____________________________________________________________________
cf291b42 91void AliFMDGainDA::Init()
92{
f631f26d 93
427e8f99 94
95 Int_t nEventsRequired = 0;
96 for(Int_t i=0;i<fEventsPerChannel.GetSize();i++)
97 {
98 Int_t nEvents = (fPulseLength.At(i)*fHighPulse) / fPulseSize.At(i);
99 fEventsPerChannel.AddAt(nEvents,i);
100 if(nEvents>nEventsRequired) nEventsRequired = nEvents * fNumberOfStripsPerChip;
101
102 }
103
54b167db 104
cf291b42 105 //8 pulser values * 128 strips * 100 samples
427e8f99 106
107
108 SetRequiredEvents(nEventsRequired);
f631f26d 109 TObjArray* detArray;
110 TObjArray* ringArray;
111 TObjArray* sectorArray;
112
113 for(UShort_t det=1;det<=3;det++) {
114 detArray = new TObjArray();
115 detArray->SetOwner();
116 fGainArray.AddAtAndExpand(detArray,det);
117 for (UShort_t ir = 0; ir < 2; ir++) {
118 Char_t ring = (ir == 0 ? 'O' : 'I');
119 UShort_t nsec = (ir == 0 ? 40 : 20);
120 UShort_t nstr = (ir == 0 ? 2 : 4);
121 ringArray = new TObjArray();
122 ringArray->SetOwner();
123 detArray->AddAtAndExpand(ringArray,ir);
124 for(UShort_t sec =0; sec < nsec; sec++) {
125 sectorArray = new TObjArray();
126 sectorArray->SetOwner();
127 ringArray->AddAtAndExpand(sectorArray,sec);
128 for(UShort_t strip = 0; strip < nstr; strip++) {
cf291b42 129 TH1S* hChannel = new TH1S(Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
130 Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
131 1024,0,1023);
f631f26d 132 hChannel->SetDirectory(0);
133 sectorArray->AddAtAndExpand(hChannel,strip);
134 }
135 }
136 }
137 }
138}
139
140//_____________________________________________________________________
141void AliFMDGainDA::AddChannelContainer(TObjArray* sectorArray,
cf291b42 142 UShort_t det ,
143 Char_t ring,
f631f26d 144 UShort_t sec,
cf291b42 145 UShort_t strip)
146{
f631f26d 147 TGraphErrors* hChannel = new TGraphErrors();
cf291b42 148 hChannel->SetName(Form("FMD%d%c[%02d,%03d]", det, ring, sec, strip));
149 hChannel->SetTitle(Form("FMD%d%c[%02d,%03d] ADC vs DAC",
150 det, ring, sec, strip));
f631f26d 151 sectorArray->AddAtAndExpand(hChannel,strip);
152}
153
154//_____________________________________________________________________
155void AliFMDGainDA::FillChannels(AliFMDDigit* digit) {
156
157 UShort_t det = digit->Detector();
158 Char_t ring = digit->Ring();
159 UShort_t sec = digit->Sector();
160 UShort_t strip = digit->Strip();
161
162
cf291b42 163 if(strip % fNumberOfStripsPerChip) return;
164 Int_t vaChip = strip / fNumberOfStripsPerChip;
165 TH1S* hChannel = GetChannelHistogram(det, ring, sec, vaChip);
f631f26d 166 hChannel->Fill(digit->Counts());
167 UpdatePulseAndADC(det,ring,sec,strip);
168}
169
170//_____________________________________________________________________
171void AliFMDGainDA::Analyse(UShort_t det,
e9c06036 172 Char_t ring,
f631f26d 173 UShort_t sec,
174 UShort_t strip) {
175 TGraphErrors* grChannel = GetChannel(det,ring,sec,strip);
176 if(!grChannel->GetN()) {
cf291b42 177 // AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
178 // det, ring , sec, strip));
f631f26d 179 return;
180 }
9c958f2b 181 TF1 fitFunc("fitFunc","pol1",-10,280);
182 fitFunc.SetParameters(100,3);
e9c06036 183 grChannel->Fit("fitFunc","Q0+","",0,fHighPulse);
f7f0b643 184
f631f26d 185 Float_t chi2ndf = 0;
9c958f2b 186 if(fitFunc.GetNDF())
187 chi2ndf = fitFunc.GetChisquare() / fitFunc.GetNDF();
f7f0b643 188
189 fOutputFile << det << ','
190 << ring << ','
191 << sec << ','
192 << strip << ','
9c958f2b 193 << fitFunc.GetParameter(1) << ','
194 << fitFunc.GetParError(1) << ','
f631f26d 195 << chi2ndf <<"\n";
196
197
3bae5d02 198 fSummaryGains.SetBinContent(fCurrentSummaryStrip,fitFunc.GetParameter(1));
199 fSummaryGains.SetBinError(fCurrentSummaryStrip,fitFunc.GetParError(1));
200 fCurrentSummaryStrip++;
f631f26d 201 if(fSaveHistograms) {
e9c06036 202 gDirectory->cd(GetSectorPath(det,ring, sec, kTRUE));
203
204 TH1F* summary = dynamic_cast<TH1F*>(gDirectory->Get("Summary"));
205 if (!summary) {
206 Int_t nStr = (ring == 'I' ? 512 : 256);
207 summary = new TH1F("Summary", Form("Summary of gains in FMD%d%c[%02d]",
208 det, ring, sec),
209 nStr, -.5, nStr-.5);
210 summary->SetXTitle("Strip");
211 summary->SetYTitle("Gain [ADC/DAC]");
212 summary->SetDirectory(gDirectory);
213 }
214 summary->SetBinContent(strip+1, fitFunc.GetParameter(1));
215 summary->SetBinError(strip+1, fitFunc.GetParError(1));
216
217 gDirectory->cd(GetStripPath(det,ring,sec,strip, kTRUE));
218 grChannel->SetName(Form("FMD%d%c[%02d,%03d]",det,ring,sec,strip));
219 // grChannel->SetDirectory(gDirectory);
220 grChannel->Write();
221 // grChannel->Write(Form("grFMD%d%c_%d_%d",det,ring,sec,strip));
222 }
f631f26d 223}
224
225//_____________________________________________________________________
3bae5d02 226void AliFMDGainDA::Terminate(TFile* diagFile)
227{
228 diagFile->cd();
229 fSummaryGains.Write();
230}
231
232//_____________________________________________________________________
e9c06036 233void AliFMDGainDA::WriteHeaderToFile()
234{
9c958f2b 235 AliFMDParameters* pars = AliFMDParameters::Instance();
236 fOutputFile.write(Form("# %s \n",pars->GetGainShuttleID()),9);
f7f0b643 237 fOutputFile.write("# Detector, "
238 "Ring, "
239 "Sector, "
e9c06036 240 "Strip, "
241 "Gain, "
242 "Error, "
f7f0b643 243 "Chi2/NDF \n",56);
f631f26d 244
245}
246
247//_____________________________________________________________________
248TH1S* AliFMDGainDA::GetChannelHistogram(UShort_t det,
e9c06036 249 Char_t ring,
f631f26d 250 UShort_t sec,
e9c06036 251 UShort_t strip)
252{
f631f26d 253
254 UShort_t Ring = 1;
255 if(ring == 'O')
256 Ring = 0;
257
258
259 TObjArray* detArray = static_cast<TObjArray*>(fGainArray.At(det));
260 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(Ring));
261 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
262 TH1S* hChannel = static_cast<TH1S*>(secArray->At(strip));
263
264 return hChannel;
265}
266
267//_____________________________________________________________________
cf291b42 268TGraphErrors* AliFMDGainDA::GetChannel(UShort_t det,
e9c06036 269 Char_t ring,
cf291b42 270 UShort_t sec,
e9c06036 271 UShort_t strip)
272{
cf291b42 273 UShort_t iring = (ring == 'O' ? 0 : 1);
274 TObjArray* detArray = static_cast<TObjArray*>(fDetectorArray.At(det));
275 TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(iring));
276 TObjArray* secArray = static_cast<TObjArray*>(ringArray->At(sec));
277 TGraphErrors* hChannel = static_cast<TGraphErrors*>(secArray->At(strip));
f631f26d 278
279 return hChannel;
280}
281
282//_____________________________________________________________________
cf291b42 283void AliFMDGainDA::UpdatePulseAndADC(UShort_t det,
284 Char_t ring,
285 UShort_t sec,
286 UShort_t strip)
287{
f7f0b643 288
427e8f99 289 AliFMDParameters* pars = AliFMDParameters::Instance();
290 UInt_t ddl, board,chip,ch;
291 pars->Detector2Hardware(det,ring,sec,strip,ddl,board,chip,ch);
292 Int_t halfring = GetHalfringIndex(det,ring,board%16);
f7f0b643 293 if(GetCurrentEvent()> (fNumberOfStripsPerChip*fEventsPerChannel.At(halfring)))
aea2699c 294 return;
cf291b42 295 if(strip % fNumberOfStripsPerChip) return;
427e8f99 296 if(((GetCurrentEvent()) % fPulseLength.At(halfring)) && GetCurrentEvent() > 0) return;
297
cf291b42 298 Int_t vaChip = strip/fNumberOfStripsPerChip;
299 TH1S* hChannel = GetChannelHistogram(det,ring,sec,vaChip);
f631f26d 300
301 if(!hChannel->GetEntries()) {
cf291b42 302 AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
303 det, ring , sec, strip));
f631f26d 304 return;
305 }
306 Double_t mean = hChannel->GetMean();
307 Double_t rms = hChannel->GetRMS();
427e8f99 308 Double_t pulse = Double_t(fCurrentPulse.At(halfring)) * fPulseSize.At(halfring);
cf291b42 309 Int_t firstBin = hChannel->GetXaxis()->GetFirst();
310 Int_t lastBin = hChannel->GetXaxis()->GetLast();
9c958f2b 311 hChannel->GetXaxis()->SetRangeUser(mean-4*rms,mean+4*rms);
f631f26d 312
cf291b42 313 mean = hChannel->GetMean();
314 rms = hChannel->GetRMS();
efd39294 315
316 hChannel->GetXaxis()->SetRange(firstBin,lastBin);
317
427e8f99 318 Int_t channelNumber = strip + (GetCurrentEvent()-1)/((fPulseLength.At(halfring)*fHighPulse)/fPulseSize.At(halfring));
f631f26d 319
320 TGraphErrors* channel = GetChannel(det,ring,sec,channelNumber);
321
427e8f99 322 channel->SetPoint(fCurrentPulse.At(halfring),pulse,mean);
323 channel->SetPointError(fCurrentPulse.At(halfring),0,rms);
f631f26d 324
325 if(fSaveHistograms) {
f7f0b643 326 gDirectory->cd(GetStripPath(det,ring,sec,channelNumber));
327 hChannel->Write(Form("%s_pulse_%03d",hChannel->GetName(),(Int_t)pulse));
328
f631f26d 329 }
330
427e8f99 331
332
333
f631f26d 334 hChannel->Reset();
335
336}
337
338//_____________________________________________________________________
cf291b42 339void AliFMDGainDA::ResetPulseAndUpdateChannel()
340{
427e8f99 341 //for(Int_t i=0; i<fCurrentPulse.GetSize();i++)
342 fCurrentPulse.Reset(0);
f631f26d 343}
344
345//_____________________________________________________________________
cf291b42 346void AliFMDGainDA::FinishEvent()
347{
427e8f99 348 for(Int_t i = 0; i<fPulseLength.GetSize();i++) {
349 if(GetCurrentEvent()>0 && (GetCurrentEvent() % fPulseLength.At(i) == 0))
350 fCurrentPulse.AddAt(fCurrentPulse.At(i)+1,i);
351
352 if(GetCurrentEvent()>0 && (GetCurrentEvent()) % fEventsPerChannel.At(i) == 0)
353 fCurrentPulse.AddAt(0,i);
354 }
f631f26d 355}
356//_____________________________________________________________________
357//
358//EOF
359//