]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDGainDA.cxx
729580af25ceda61eae54c2e79459d54d37dcbf6
[u/mrichter/AliRoot.git] / FMD / AliFMDGainDA.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    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 */
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 //
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 #include "AliFMDAltroMapping.h"
41
42 //_____________________________________________________________________
43 ClassImp(AliFMDGainDA)
44 #if 0 // Do not delete - here to let Emacs indent properly
45 ;
46 #endif
47
48 //_____________________________________________________________________
49 AliFMDGainDA::AliFMDGainDA() 
50   : AliFMDBaseDA(),
51     fGainArray(),
52     fHighPulse(256), 
53     fEventsPerChannel(10),
54     fCurrentPulse(10),
55     fCurrentChannel(10),
56     fNumberOfStripsPerChip(128),
57     fSummaryGains("GainsSummary","Summary of gains",51200,0,51200),
58     fCurrentSummaryStrip(1)
59 {
60   fCurrentPulse.Reset(0);
61   fCurrentChannel.Reset(0);
62   fOutputFile.open("gains.csv");
63   fGainArray.SetOwner(); 
64 }
65
66 //_____________________________________________________________________
67 AliFMDGainDA::AliFMDGainDA(const AliFMDGainDA & gainDA) 
68   :  AliFMDBaseDA(gainDA),
69      fGainArray(gainDA.fGainArray),
70      fHighPulse(gainDA.fHighPulse),
71      fEventsPerChannel(gainDA.fEventsPerChannel),
72      fCurrentPulse(gainDA.fCurrentPulse),
73      fCurrentChannel(gainDA.fCurrentChannel),
74      fNumberOfStripsPerChip(gainDA.fNumberOfStripsPerChip),
75      fSummaryGains(gainDA.fSummaryGains),
76      fCurrentSummaryStrip(gainDA.fCurrentSummaryStrip)
77 {  
78   fCurrentPulse.Reset(0);
79   fCurrentChannel.Reset(0);
80 }
81
82 //_____________________________________________________________________
83 AliFMDGainDA::~AliFMDGainDA() 
84 {
85 }
86
87 //_____________________________________________________________________
88 void AliFMDGainDA::Init() 
89 {
90   
91  
92   Int_t nEventsRequired = 0;
93   
94   //for(UShort_t det=1; det<=3;det++) {
95   //  UShort_t firstring = (det == 1 ? 1 : 0);
96   //  for(UShort_t iring = firstring; iring <=1;iring++) {
97   //    Char_t ring = (iring == 1 ? 'I' : 'O');
98   //    for(UShort_t board =0 ; board <=1; board++) {
99   //    Int_t idx = GetHalfringIndex(det,ring,board);
100   for(Int_t idx = 0;idx<fEventsPerChannel.GetSize();idx++)
101     {
102       
103       Int_t nEvents = 0;
104       if(fPulseSize.At(idx))
105         nEvents = (fPulseLength.At(idx)*fHighPulse) / fPulseSize.At(idx);
106       fEventsPerChannel.AddAt(nEvents,idx);
107       if(nEvents>nEventsRequired) nEventsRequired = nEvents * fNumberOfStripsPerChip;
108       
109     }
110   //}
111   // }
112   
113   //8 pulser values * 128 strips * 100 samples
114   
115   
116   SetRequiredEvents(nEventsRequired); 
117   
118   TObjArray* detArray;
119   TObjArray* ringArray;
120   TObjArray* sectorArray;
121   
122   for(UShort_t det=1;det<=3;det++) {
123     detArray = new TObjArray();
124     detArray->SetOwner();
125     fGainArray.AddAtAndExpand(detArray,det);
126     for (UShort_t ir = 0; ir < 2; ir++) {
127       Char_t   ring = (ir == 0 ? 'O' : 'I');
128       UShort_t nsec = (ir == 0 ? 40  : 20);
129       UShort_t nstr = (ir == 0 ? 2 : 4);
130       ringArray = new TObjArray();
131       ringArray->SetOwner();
132       detArray->AddAtAndExpand(ringArray,ir);
133       for(UShort_t sec =0; sec < nsec;  sec++)  {
134         sectorArray = new TObjArray();
135         sectorArray->SetOwner();
136         ringArray->AddAtAndExpand(sectorArray,sec);
137         for(UShort_t strip = 0; strip < nstr; strip++) {
138           TH1S* hChannel = new TH1S(Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
139                                     Form("hFMD%d%c_%d_%d",det,ring,sec,strip),
140                                     1024,0,1023);
141           hChannel->SetDirectory(0);
142           sectorArray->AddAtAndExpand(hChannel,strip);
143         }
144       }
145     }
146   }
147 }
148
149 //_____________________________________________________________________
150 void AliFMDGainDA::AddChannelContainer(TObjArray* sectorArray, 
151                                        UShort_t det  , 
152                                        Char_t   ring,  
153                                        UShort_t sec, 
154                                        UShort_t strip) 
155 {  
156   TGraphErrors* hChannel  = new TGraphErrors();
157   hChannel->SetName(Form("FMD%d%c[%02d,%03d]", det, ring, sec, strip));
158   hChannel->SetTitle(Form("FMD%d%c[%02d,%03d] ADC vs DAC", 
159                           det, ring, sec, strip));
160   sectorArray->AddAtAndExpand(hChannel,strip);
161 }
162
163 //_____________________________________________________________________
164 void AliFMDGainDA::FillChannels(AliFMDDigit* digit) {
165
166   UShort_t det   = digit->Detector();
167   Char_t   ring  = digit->Ring();
168   UShort_t sec   = digit->Sector();
169   UShort_t strip = digit->Strip();
170   
171   //Strip is always seen as the first in a VA chip. All other strips are junk
172   if(strip % fNumberOfStripsPerChip) return;
173   Int_t vaChip   = strip / fNumberOfStripsPerChip; 
174   TH1S* hChannel = GetChannelHistogram(det, ring, sec, vaChip);
175   
176   hChannel->Fill(digit->Counts());
177   UpdatePulseAndADC(det,ring,sec,strip);
178 }
179
180 //_____________________________________________________________________
181 void AliFMDGainDA::Analyse(UShort_t det, 
182                            Char_t   ring, 
183                            UShort_t sec, 
184                            UShort_t strip) {
185   TGraphErrors* grChannel = GetChannel(det,ring,sec,strip);
186   if(!grChannel->GetN()) {
187     // AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
188     //                 det, ring , sec, strip));
189     return;
190   }
191   TF1 fitFunc("fitFunc","pol1",-10,280); 
192   fitFunc.SetParameters(100,3);
193   grChannel->Fit("fitFunc","Q0+","",0,fHighPulse);
194      
195   Float_t gain    = -1;
196   Float_t error   = -1; 
197   Float_t chi2ndf = -1;
198   if((fitFunc.GetParameter(1)) == (fitFunc.GetParameter(1))) {
199     gain    = fitFunc.GetParameter(1);
200     error   = fitFunc.GetParError(1);
201     if(fitFunc.GetNDF())
202       chi2ndf = fitFunc.GetChisquare() / fitFunc.GetNDF();
203   }
204   
205   fOutputFile << det                         << ','
206               << ring                        << ','
207               << sec                         << ','
208               << strip                       << ','
209               << gain                        << ','
210               << error                       << ','
211               << chi2ndf                     <<"\n";
212   
213   //due to RCU trouble, first strips on VAs are excluded
214   if(strip%128 != 0) {
215     
216     fSummaryGains.SetBinContent(fCurrentSummaryStrip,fitFunc.GetParameter(1));
217     fSummaryGains.SetBinError(fCurrentSummaryStrip,fitFunc.GetParError(1));
218     
219     fCurrentSummaryStrip++;
220   }
221   if(fSaveHistograms) {
222     gDirectory->cd(GetSectorPath(det,ring, sec, kTRUE));
223     
224     TH1F* summary = dynamic_cast<TH1F*>(gDirectory->Get("Summary"));
225     if (!summary) { 
226       Int_t nStr = (ring == 'I' ? 512 : 256);
227       summary = new TH1F("Summary", Form("Summary of gains in FMD%d%c[%02d]", 
228                                          det, ring, sec), 
229                          nStr, -.5, nStr-.5);
230       summary->SetXTitle("Strip");
231       summary->SetYTitle("Gain [ADC/DAC]");
232       summary->SetDirectory(gDirectory);
233     }
234     summary->SetBinContent(strip+1, fitFunc.GetParameter(1));
235     summary->SetBinError(strip+1, fitFunc.GetParError(1));
236     
237     gDirectory->cd(GetStripPath(det,ring,sec,strip, kTRUE));
238     grChannel->SetName(Form("FMD%d%c[%02d,%03d]",det,ring,sec,strip));
239     // grChannel->SetDirectory(gDirectory);
240     grChannel->Write();
241     // grChannel->Write(Form("grFMD%d%c_%d_%d",det,ring,sec,strip));
242   }  
243 }
244
245 //_____________________________________________________________________
246 void AliFMDGainDA::Terminate(TFile* diagFile)
247 {
248   if(diagFile) {
249     diagFile->cd();
250     fSummaryGains.Write();
251   }
252 }
253
254 //_____________________________________________________________________
255 void AliFMDGainDA::WriteHeaderToFile() 
256 {
257   AliFMDParameters* pars       = AliFMDParameters::Instance();
258   fOutputFile.write(Form("# %s \n",pars->GetGainShuttleID()),9);
259   fOutputFile.write("# Detector, "
260                     "Ring, "
261                     "Sector, "
262                     "Strip, "
263                     "Gain, "
264                     "Error, "
265                     "Chi2/NDF \n",56);
266   
267 }
268
269 //_____________________________________________________________________
270 TH1S* AliFMDGainDA::GetChannelHistogram(UShort_t det, 
271                                         Char_t   ring, 
272                                         UShort_t sec, 
273                                         UShort_t strip) 
274 {
275   
276   UShort_t  Ring = 1;
277   if(ring == 'O')
278     Ring = 0;
279   
280   
281   TObjArray* detArray  = static_cast<TObjArray*>(fGainArray.At(det));
282   TObjArray* ringArray = static_cast<TObjArray*>(detArray->At(Ring));
283   TObjArray* secArray  = static_cast<TObjArray*>(ringArray->At(sec));
284   TH1S* hChannel       = static_cast<TH1S*>(secArray->At(strip));
285   
286   return hChannel;
287 }
288
289 //_____________________________________________________________________
290 TGraphErrors* AliFMDGainDA::GetChannel(UShort_t det, 
291                                        Char_t   ring, 
292                                        UShort_t sec, 
293                                        UShort_t strip) 
294 {  
295   UShort_t      iring     = (ring == 'O' ? 0 : 1);
296   TObjArray*    detArray  = static_cast<TObjArray*>(fDetectorArray.At(det));
297   TObjArray*    ringArray = static_cast<TObjArray*>(detArray->At(iring));
298   TObjArray*    secArray  = static_cast<TObjArray*>(ringArray->At(sec));
299   TGraphErrors* hChannel  = static_cast<TGraphErrors*>(secArray->At(strip));
300   
301   return hChannel;
302 }
303
304 //_____________________________________________________________________
305 void AliFMDGainDA::UpdatePulseAndADC(UShort_t det, 
306                                      Char_t ring, 
307                                      UShort_t sec, 
308                                      UShort_t strip) 
309 {
310   
311   AliFMDParameters* pars = AliFMDParameters::Instance();
312   // UInt_t ddl, board,chip,ch;
313   UShort_t board = pars->GetAltroMap()->Sector2Board(ring, sec);
314   // pars->Detector2Hardware(det,ring,sec,strip,ddl,board,chip,ch);
315   /// pars->GetAltroMap()->Strip2Channel(
316   Int_t halfring = GetHalfringIndex(det,ring,board/16);
317   
318   if(GetCurrentEvent()> (fNumberOfStripsPerChip*fEventsPerChannel.At(halfring)))
319     return;
320   if(strip % fNumberOfStripsPerChip) return;
321   if(((GetCurrentEvent()) % fPulseLength.At(halfring)) 
322      && GetCurrentEvent() > 0) return;
323      
324   Int_t vaChip = strip/fNumberOfStripsPerChip; 
325   TH1S* hChannel = GetChannelHistogram(det,ring,sec,vaChip);
326   
327   if(!hChannel->GetEntries()) {
328     AliWarning(Form("No entries for FMD%d%c, sector %d, strip %d",
329                     det, ring , sec, strip));
330     return;
331   }
332   Double_t mean      = hChannel->GetMean();
333   Double_t rms       = hChannel->GetRMS();
334   Double_t pulse     = (Double_t(fCurrentPulse.At(halfring)) 
335                         * fPulseSize.At(halfring));
336   Int_t    firstBin  = hChannel->GetXaxis()->GetFirst();
337   Int_t    lastBin   = hChannel->GetXaxis()->GetLast();
338   hChannel->GetXaxis()->SetRangeUser(mean-4*rms,mean+4*rms);
339   
340   mean               = hChannel->GetMean();
341   rms                = hChannel->GetRMS();
342   
343   hChannel->GetXaxis()->SetRange(firstBin,lastBin);
344   
345   Int_t channelNumber = (strip + 
346                          (GetCurrentEvent()-1)
347                          / ((fPulseLength.At(halfring)*fHighPulse)
348                             / fPulseSize.At(halfring))); 
349   
350   TGraphErrors* channel = GetChannel(det,ring,sec,channelNumber);
351   
352   channel->SetPoint(fCurrentPulse.At(halfring),pulse,mean);
353   channel->SetPointError(fCurrentPulse.At(halfring),0,rms);
354   
355   if(fSaveHistograms) {
356     gDirectory->cd(GetStripPath(det,ring,sec,channelNumber));
357     hChannel->Write(Form("%s_pulse_%03d",hChannel->GetName(),(Int_t)pulse));
358     
359   }
360     
361   hChannel->Reset();
362   
363 }
364
365 //_____________________________________________________________________
366 void AliFMDGainDA::ResetPulseAndUpdateChannel() 
367 {  
368   fCurrentPulse.Reset(0); 
369 }
370
371 //_____________________________________________________________________
372 void AliFMDGainDA::FinishEvent() 
373 {
374   for(UShort_t det=1; det<=3;det++) {
375     UShort_t firstring = (det == 1 ? 1 : 0);
376     for(UShort_t iring = firstring; iring <=1;iring++) {
377       Char_t ring = (iring == 1 ? 'I' : 'O');
378       for(UShort_t board =0 ; board <=1; board++) {
379         Int_t idx = GetHalfringIndex(det,ring,board);
380         
381         if( !fPulseLength.At(idx) || !fEventsPerChannel.At(idx))
382           continue;
383         if(GetCurrentEvent()>0 && ((GetCurrentEvent() % fPulseLength.At(idx)) == 0))
384           fCurrentPulse.AddAt(fCurrentPulse.At(idx)+1,idx);
385         
386         if(GetCurrentEvent()>0 && ((GetCurrentEvent()) % fEventsPerChannel.At(idx)) == 0)
387           fCurrentPulse.AddAt(0,idx);
388       }
389     }
390   }
391 }
392 //_____________________________________________________________________
393 //
394 //EOF
395 //