]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDPedestalDA.cxx
Always export to FES the regional configuration file
[u/mrichter/AliRoot.git] / FMD / AliFMDPedestalDA.cxx
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 //
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.
28 //
29
30 #include "AliFMDPedestalDA.h"
31 #include "AliFMDAltroMapping.h"
32 #include "iostream"
33 #include "fstream"
34 #include "AliLog.h"
35 #include "TF1.h"
36 #include "TObject.h"
37 #include "TMath.h"
38
39 //_____________________________________________________________________
40 ClassImp(AliFMDPedestalDA)
41
42 //_____________________________________________________________________
43 AliFMDPedestalDA::AliFMDPedestalDA() : AliFMDBaseDA(),
44   fCurrentChannel(1),
45   fPedSummary("PedestalSummary","pedestals",51200,0,51200),
46   fNoiseSummary("NoiseSummary","noise",51200,0,51200),
47   fZSfileFMD1(),
48   fZSfileFMD2(),
49   fZSfileFMD3(), 
50   fMinTimebin(3 * 4 * 3 * 16), // 3 ddls, 4 FECs, 3 Altros, 16 channels
51   fMaxTimebin(3 * 4 * 3 * 16)  // 3 ddls, 4 FECs, 3 Altros, 16 channels
52 {
53   // Default constructor 
54   fOutputFile.open("peds.csv");
55   fZSfileFMD1.open("ddl3072.csv");
56   fZSfileFMD2.open("ddl3073.csv");
57   fZSfileFMD3.open("ddl3074.csv");  
58 }
59
60 //_____________________________________________________________________
61 AliFMDPedestalDA::AliFMDPedestalDA(const AliFMDPedestalDA & pedDA) : 
62   AliFMDBaseDA(pedDA),
63   fCurrentChannel(1),
64   fPedSummary("PedestalSummary","pedestals",51200,0,51200),
65   fNoiseSummary("NoiseSummary","noise",51200,0,51200),
66   fZSfileFMD1(),
67   fZSfileFMD2(),
68   fZSfileFMD3(),
69   fMinTimebin(pedDA.fMinTimebin),
70   fMaxTimebin(pedDA.fMaxTimebin)
71 {
72   // Copy constructor 
73 }
74
75 //_____________________________________________________________________
76 AliFMDPedestalDA::~AliFMDPedestalDA() 
77
78   // Destructor.
79 }
80
81 //_____________________________________________________________________
82 void AliFMDPedestalDA::Init() 
83
84   // Initialise 
85   SetRequiredEvents(1000);
86   fMinTimebin.Reset(1024);
87   fMaxTimebin.Reset(-1);
88 }
89
90 //_____________________________________________________________________
91 void AliFMDPedestalDA::AddChannelContainer(TObjArray* sectorArray, 
92                                            UShort_t det, 
93                                            Char_t   ring, 
94                                            UShort_t sec, 
95                                            UShort_t strip) 
96 {
97   // Add a channel to the containers.
98   //
99   // Parameters:
100   //     sectorArray  Array of sectors
101   //     det          Detector 
102   //     ring         Ring
103   //     sec          Sector 
104   //     strip        Strip
105   AliFMDParameters* pars        = AliFMDParameters::Instance();
106   UInt_t            samples     = pars->GetSampleRate(det, ring, sec, strip);
107   TObjArray*        sampleArray = new TObjArray(samples);
108   sampleArray->SetOwner();
109   for (UInt_t sample = 0; sample < samples; sample++) {
110     TH1S* hSample = new TH1S(Form("FMD%d%c[%02d,03%d]_%d",
111                                   det,ring,sec,strip,sample),
112                              Form("FMD%d%c[%02d,%03%d]_%d",
113                                   det,ring,sec,strip),
114                              1024,-.5,1023.5);
115     hSample->SetXTitle("ADC");
116     hSample->SetYTitle("Events");
117     hSample->SetDirectory(0);
118     sampleArray->AddAt(hSample, sample);
119   }
120   sectorArray->AddAtAndExpand(sampleArray, strip);
121 }
122
123 //_____________________________________________________________________
124 void AliFMDPedestalDA::FillChannels(AliFMDDigit* digit) 
125 {
126   // Fill ADC values from a digit into the corresponding histogram.
127   //
128   // Parameters:
129   //    digit Digit to fill ADC values for.
130   UShort_t det   = digit->Detector();
131   Char_t   ring  = digit->Ring();
132   UShort_t sec   = digit->Sector();
133   UShort_t strip = digit->Strip();
134   
135   AliFMDParameters* pars     = AliFMDParameters::Instance();
136   UInt_t            samples  = pars->GetSampleRate(det, ring, sec, strip);
137   for (UInt_t sample = 0; sample < samples; sample++) {
138     TH1S* hSample = GetChannel(det, ring, sec, strip, sample);
139     hSample->Fill(digit->Count(sample));
140   }
141   
142 }
143
144 //_____________________________________________________________________
145 void AliFMDPedestalDA::Analyse(UShort_t det, 
146                                Char_t   ring, 
147                                UShort_t sec, 
148                                UShort_t strip) 
149 {
150   // Analyse a strip.  That is, compute the mean and spread of the ADC
151   // spectra for all strips.  Also output on files the values. 
152   //
153   // Parameters:
154   //     det   Detector
155   //     ring  Ring
156   //     sec   Sector 
157   //     strip Strip.
158   AliFMDParameters* pars     = AliFMDParameters::Instance();
159   // Float_t           factor   = pars->GetPedestalFactor();
160   UInt_t            samples  = pars->GetSampleRate(det, ring, sec, strip);
161   for (UShort_t sample = 0; sample < samples; sample++) {
162   
163     TH1S* hChannel = GetChannel(det, ring, sec, strip,sample);
164     if(hChannel->GetEntries() == 0) {
165       //AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
166       //            det,ring,sec,strip));
167       return;
168     }
169     
170     AliDebug(50, Form("Fitting FMD%d%c_%d_%d with %d entries",
171                       det,ring,sec,strip, hChannel->GetEntries()));
172     TF1 fitFunc("fitFunc","gausn",0,300);
173     fitFunc.SetParameters(100,100,1);
174     hChannel->Fit("fitFunc","Q0+","",10,200);
175     
176     Float_t mean = hChannel->GetMean();
177     Float_t rms  = hChannel->GetRMS();
178     
179     
180     
181     hChannel->GetXaxis()->SetRangeUser(mean-5*rms,mean+5*rms);
182     
183     mean = hChannel->GetMean();
184     rms  = hChannel->GetRMS();
185   
186     
187     UShort_t ddl, board, altro, channel;
188     UShort_t timebin;
189     
190     pars->Detector2Hardware(det,ring,sec,strip,sample,
191                             ddl,board,altro,channel,timebin);
192     Int_t idx = HWIndex(ddl, board, altro, channel);
193     if (idx >= 0) { 
194       fMinTimebin[idx] = TMath::Min(Short_t(timebin),   fMinTimebin[idx]);
195       fMaxTimebin[idx] = TMath::Max(Short_t(timebin+1), fMaxTimebin[idx]);
196     }
197     
198     std::ostream* zsFile = 0;
199     switch(det) {
200     case 1:  zsFile = &fZSfileFMD1; break;
201     case 2:  zsFile = &fZSfileFMD2; break;
202     case 3:  zsFile = &fZSfileFMD3; break;
203     default: AliWarning("Unknown sample!"); break;
204       
205     }
206     *zsFile << board   << ',' 
207             << altro   << ',' 
208             << channel << ',' 
209             << timebin << ','
210             << mean    << ',' 
211             << rms     << "\n"; break;
212     
213     Float_t chi2ndf = 0;
214     
215     
216     if(fitFunc.GetNDF())
217       chi2ndf = fitFunc.GetChisquare() / fitFunc.GetNDF();
218     
219     
220     Int_t sampleToWrite = 2;
221     if      (samples == 2) sampleToWrite = 1;
222     else if (samples <  2) sampleToWrite = 0;
223     
224     if(sample != sampleToWrite) continue;
225     
226     
227     fOutputFile << det                         << ','
228                 << ring                        << ','
229                 << sec                         << ','
230                 << strip                       << ','
231                 << mean                        << ','
232                 << rms                         << ','
233                 << fitFunc.GetParameter(1)     << ','
234                 << fitFunc.GetParameter(2)     << ','
235                 << chi2ndf                     <<"\n";
236     
237     if(fSaveHistograms  ) {
238       gDirectory->cd(GetSectorPath(det, ring, sec, kTRUE));
239       TH1F* sumPed   = dynamic_cast<TH1F*>(gDirectory->Get("Pedestals"));
240       TH1F* sumNoise = dynamic_cast<TH1F*>(gDirectory->Get("Noise"));
241       Int_t nStr = (ring == 'I' ? 512 : 256);
242       if (!sumPed) {
243         sumPed = new TH1F("Pedestals", 
244                           Form("Summary of pedestals in FMD%d%c[%02d]", 
245                                det, ring, sec), 
246                           nStr, -.5, nStr-.5);
247         sumPed->SetXTitle("Strip");
248         sumPed->SetYTitle("Pedestal [ADC]");
249         sumPed->SetDirectory(gDirectory);
250       }
251       if (!sumNoise) { 
252         sumNoise = new TH1F("Noise", 
253                             Form("Summary of noise in FMD%d%c[%02d]", 
254                                  det, ring, sec), 
255                             nStr, -.5, nStr-.5);
256         sumNoise->SetXTitle("Strip");
257         sumNoise->SetYTitle("Noise [ADC]");
258         
259         sumNoise->SetDirectory(gDirectory);
260       }
261       sumPed->SetBinContent(strip+1, mean);
262       sumPed->SetBinError(strip+1, rms);
263       sumNoise->SetBinContent(strip+1, rms);
264       
265       if(sumNoise->GetEntries() == nStr)
266         sumNoise->Write(sumNoise->GetName(),TObject::kOverwrite);
267       if(sumPed->GetEntries() == nStr)
268         sumPed->Write(sumPed->GetName(),TObject::kOverwrite);
269       
270       fPedSummary.SetBinContent(fCurrentChannel,mean);
271       
272       fNoiseSummary.SetBinContent(fCurrentChannel,rms);
273       fCurrentChannel++;
274       
275       gDirectory->cd(GetStripPath(det, ring, sec, strip, kTRUE));
276       hChannel->GetXaxis()->SetRange(1,1024);
277       
278       hChannel->Write();
279     }
280   }
281 }
282
283 //_____________________________________________________________________
284 void AliFMDPedestalDA::Terminate(TFile* diagFile) 
285 {
286   // Called at the end of a job.  Fills in missing time-bins and
287   // closes output files
288   if(fSaveHistograms) {
289     diagFile->cd();
290     
291     fPedSummary.Write();
292     fNoiseSummary.Write();
293   }
294   AliFMDAltroMapping* map = AliFMDParameters::Instance()->GetAltroMap();
295   if(fZSfileFMD1.is_open()) FillinTimebins(fZSfileFMD1, map->Detector2DDL(1));
296   if(fZSfileFMD2.is_open()) FillinTimebins(fZSfileFMD2, map->Detector2DDL(2));
297   if(fZSfileFMD3.is_open()) FillinTimebins(fZSfileFMD3, map->Detector2DDL(3));
298   
299 }
300
301 //_____________________________________________________________________
302 void AliFMDPedestalDA::FillinTimebins(std::ofstream& out, UShort_t /*ddl*/)
303 {
304 #if 0
305   unsigned short  boards[] = { 0x0, 0x1, 0x10, 0x11, 0xFFFF };
306   unsigned short* board    = boards;
307   while ((*boards) != 0xFFFF) { 
308     for (UShort_t altro = 0; altro < 3; altro++) { 
309       for (UShort_t channel = 0; channel < 16; channel++) { 
310         Int_t idx = HWIndex(ddl, *board, altro, channel);
311         if (idx < 0) { 
312           AliWarning(Form("Invalid index for %4d/0x%02x/0x%x/0x%x: %d", 
313                           ddl, *board, altro, channel, idx));
314           continue;
315         }
316         Short_t min = fMinTimebin[idx];
317         Short_t max = fMaxTimebin[idx];
318         
319         // Channel not seen at all. 
320         if (min > 1023 || max < 0) continue;
321
322         out << "# Extra timebins for 0x" << std::hex 
323             << board << ',' << altro << ',' << channel 
324             << " got time-bins " << min << " to " << max-1 
325             << std::dec << std::endl;
326         
327         for (UShort_t t = 15; t < min; t++) 
328           // Write a phony line 
329           out << board << "," << altro << "," << channel << "," 
330               << t << "," << 1023 << "," << 0 << std::endl;
331
332         for (UShort_t t = max; t < 1024; t++) 
333           // Write a phony line 
334           out << board << "," << altro << "," << channel << "," 
335               << t << "," << 1023 << "," << 0 << std::endl;
336       } // channel loop
337     } // altro loop
338   } // board loop
339   // Write trailer, and close 
340 #endif
341   out.write("# EOF\n", 6);
342   out.close();
343 }
344
345 //_____________________________________________________________________
346 void AliFMDPedestalDA::WriteHeaderToFile() 
347 {
348   // Write headers to output files 
349   AliFMDParameters* pars       = AliFMDParameters::Instance();
350   fOutputFile.write(Form("# %s \n",pars->GetPedestalShuttleID()),13);
351   fOutputFile.write("# Detector, "
352                     "Ring, "
353                     "Sector, "
354                     "Strip, "
355                     "Pedestal, "
356                     "Noise, "
357                     "Mu, "
358                     "Sigma, "
359                     "Chi2/NDF \n", 71);
360
361   std::ostream* zss[] = { &fZSfileFMD1, &fZSfileFMD2, &fZSfileFMD3, 0 };
362   for (size_t i = 0; i < 3; i++) 
363     *(zss[i]) << "# FMD 1 pedestals \n"
364               << "# board, "
365               << "altro, "
366               << "channel, "
367               << "timebin, "
368               << "pedestal, "
369               << "noise\n";
370 }
371
372 //_____________________________________________________________________
373 TH1S* AliFMDPedestalDA::GetChannel(UShort_t det, 
374                                    Char_t   ring, 
375                                    UShort_t sec, 
376                                    UShort_t strip,
377                                    UInt_t   sample) 
378 {
379   // Get the histogram corresponding to a strip sample.
380   // 
381   // Parameters:
382   //     det    Detector
383   //     ring   Ring
384   //     sec    Sector
385   //     strip  Strip
386   //     sample Sample
387   //
388   // Return:
389   //     ADC spectra of a strip.
390   UShort_t   iring       = (ring == 'O' ? 0 : 1);
391   TObjArray* detArray    = static_cast<TObjArray*>(fDetectorArray.At(det));
392   TObjArray* ringArray   = static_cast<TObjArray*>(detArray->At(iring));
393   TObjArray* secArray    = static_cast<TObjArray*>(ringArray->At(sec));
394   TObjArray* sampleArray = static_cast<TObjArray*>(secArray->At(strip));
395   TH1S*      hSample     = static_cast<TH1S*>(sampleArray->At(sample));
396   return hSample;
397   
398 }
399
400 //_____________________________________________________________________
401 //
402 //EOF
403 //