]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDCalib.cxx
High value of mean and sigma in case of no entries in the pad sigma calculations
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDCalib.cxx
1 //  **************************************************************************
2 //  * Copyright(c) 1998-1999, 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 #include "AliHMPIDCalib.h" //class header
17 #include "AliHMPIDParam.h" //class header
18 #include "AliHMPIDRawStream.h" //class header
19 #include <fstream>
20 #include <TTree.h>
21
22
23
24 ClassImp(AliHMPIDCalib) 
25
26
27 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28 AliHMPIDCalib::AliHMPIDCalib():
29 faddl(0x0),
30 fsq(0x0),
31 fsq2(0x0),
32 fnpc(0x0),
33 fpedQ0(0x0),
34 fErr(0x0),
35 fPadAdc(0x0),
36 fIsPad(0x0),
37 fFile(0x0),
38 fLdcId(0),
39 fTimeStamp(0),
40 fRunNum(0),
41 fSigCut(0),
42 fWritePads(0),
43 fnDDLInStream(0x0),
44 fnDDLOutStream(0x0),
45 fLargeHisto(kFALSE),
46 fSelectDDL(0)  
47 {
48   //
49   //constructor
50   //
51   faddl = new Bool_t[AliHMPIDRawStream::kNDDL];
52   Int_t nPads =  (AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPcx+1)*(AliHMPIDParam::kMaxPcy+1);
53  
54   fpedQ0 = new Int_t***[AliHMPIDRawStream::kNDDL+1];
55   fsq2   = new Float_t ***[AliHMPIDRawStream::kNDDL+1];
56   fsq    = new Float_t ***[AliHMPIDRawStream::kNDDL+1];
57   fnpc   = new Int_t ***[AliHMPIDRawStream::kNDDL+1];
58     fErr = new Int_t*[AliHMPIDRawStream::kNDDL+1];
59    
60   fnDDLInStream  = new Int_t[AliHMPIDRawStream::kNDDL+1];
61   fnDDLOutStream = new Int_t[AliHMPIDRawStream::kNDDL+1];
62
63   
64   for(Int_t iDDL=0;iDDL<AliHMPIDRawStream::kNDDL+1;iDDL++) {
65     
66       fErr[iDDL] = new Int_t[AliHMPIDRawStream::kSumErr+1];
67     fpedQ0[iDDL] = new Int_t**[AliHMPIDRawStream::kNRows+1];
68        fsq[iDDL] = new Float_t**[AliHMPIDRawStream::kNRows+1];
69       fsq2[iDDL] = new Float_t**[AliHMPIDRawStream::kNRows+1];
70       fnpc[iDDL] = new Int_t**[AliHMPIDRawStream::kNRows+1];
71       
72       for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++)  {
73       
74        fpedQ0[iDDL][iRow] = new Int_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
75           fsq[iDDL][iRow] = new Float_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
76          fsq2[iDDL][iRow] = new Float_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
77          fnpc[iDDL][iRow] = new Int_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
78       
79         for(Int_t iDil=1;iDil<AliHMPIDRawStream::kNDILOGICAdd+1;iDil++){
80       
81          fpedQ0[iDDL][iRow][iDil] = new Int_t[AliHMPIDRawStream::kNPadAdd+1];
82            fsq2[iDDL][iRow][iDil] = new Float_t[AliHMPIDRawStream::kNPadAdd+1];
83             fsq[iDDL][iRow][iDil] = new Float_t[AliHMPIDRawStream::kNPadAdd+1];
84            fnpc[iDDL][iRow][iDil] = new Int_t[AliHMPIDRawStream::kNPadAdd+1];
85           }//iDil
86       }//iRow
87    }//iDDL
88     
89    for(Int_t iDDL=0;iDDL<AliHMPIDRawStream::kNDDL+1;iDDL++) {
90         
91      fnDDLInStream[iDDL]=-1;
92      fnDDLOutStream[iDDL]=-1;
93       
94      for(Int_t iErr=0;iErr<AliHMPIDRawStream::kSumErr+1;iErr++)  {fErr[iDDL][iErr]=0;}
95          
96      for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++) {
97         for(Int_t iDil=1;iDil<AliHMPIDRawStream::kNDILOGICAdd+1;iDil++) {
98           for(Int_t iPad=1;iPad<AliHMPIDRawStream::kNPadAdd+1;iPad++) {
99             fpedQ0[iDDL][iRow][iDil][iPad]=0;
100                fsq[iDDL][iRow][iDil][iPad]=0;
101               fsq2[iDDL][iRow][iDil][iPad]=0;
102               fnpc[iDDL][iRow][iDil][iPad]=0;
103         }//iPad
104       }//iDil
105      }//iRow
106    }//iDDL
107     
108   fPadAdc=new TH1I*[nPads];  
109   fIsPad=new Bool_t[nPads];  
110   for(Int_t np=0;np<nPads;np++) {fPadAdc[np]=0x0;   fIsPad[np]=kFALSE;}
111   fWritePads=kFALSE;
112
113
114   Init();
115 }
116 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
117 AliHMPIDCalib::~AliHMPIDCalib()
118 {
119   //
120   //destructor
121   //
122   if (faddl)     { delete [] faddl;   faddl = 0x0;  } 
123   if (fPadAdc)   { delete [] fPadAdc; fPadAdc=0x0;  }  
124   if (fIsPad)    { delete [] fIsPad;  fIsPad=0x0;   }  
125   if (fFile)     { delete    fFile;   fFile=0x0;    }  
126  
127   for(Int_t iErr=0;iErr<AliHMPIDRawStream::kSumErr+1;iErr++) { delete [] fErr[iErr];}  delete [] fErr;
128   
129   for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
130    for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++)         
131      for(Int_t iDil=1;iDil<AliHMPIDRawStream::kNDILOGICAdd+1;iDil++)
132       {
133          delete [] fpedQ0[iDDL][iRow][iDil]; //del iPad
134          delete []    fsq[iDDL][iRow][iDil]; //del iPad
135          delete []   fsq2[iDDL][iRow][iDil]; //del iPad
136          delete []   fnpc[iDDL][iRow][iDil]; //del iPad
137        }
138    for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
139      for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++)         
140       {
141         delete [] fpedQ0[iDDL][iRow];  //del iRow
142           delete []  fsq[iDDL][iRow];  //del iRow
143           delete [] fsq2[iDDL][iRow];  //del iRow
144           delete [] fnpc[iDDL][iRow];  //del iRow
145         }
146        
147    for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
148    {   
149        delete [] fpedQ0[iDDL];        //del iRow
150          delete [] fsq2[iDDL];        //del iRow
151          delete []  fsq[iDDL];        //del iRow
152          delete [] fnpc[iDDL];        //del iRow
153      }
154        
155    delete [] fpedQ0;
156    delete [] fsq2;
157    delete [] fsq;
158    delete [] fnpc;
159     
160   fpedQ0=0;    
161     fsq2=0;
162      fsq=0;
163     fnpc=0;
164     
165   fLdcId=0;
166   fTimeStamp=0;
167   fRunNum=0;
168   fSigCut=0;
169   fWritePads=0;
170   fLargeHisto=kFALSE;
171   fSelectDDL=0;
172 }//dtor
173 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
174 void AliHMPIDCalib::Init()
175 {
176   //
177   //Init the q calc.
178   //Arguments: none
179   //Return: none
180   //
181     fSigCut=3;  
182     for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
183       {
184          for(Int_t ierr=0; ierr <AliHMPIDRawStream::kSumErr ; ierr++) {
185             fErr[iDDL][ierr]=0;
186             }
187         
188         faddl[iDDL]=kFALSE;
189       }//DDL
190 }//Init()
191 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
192 void AliHMPIDCalib::SetRunParams(ULong_t runNum,Int_t timeStamp, Int_t ldcId)
193 {
194   //  
195   //Set run parameters for the Pedestal and Error Files
196   //Arguments: run number, time stamp and LDC Id
197   //Returns: none
198   //
199   fRunNum=(Int_t)runNum;
200   fTimeStamp=timeStamp;
201   fLdcId=ldcId;
202 }//SetRunParams()
203 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
204 void AliHMPIDCalib::SetSigCutFromFile(Char_t* name)
205 {
206   //
207   //Set Sigma Cut from the file on the LDC, if the input file is not present default value is set!
208   //Arguments: the name of the SigmaCut file on the LDC
209   //Returns: none
210   //
211   Int_t nSigCut=0;
212   ifstream infile(name);
213   if(!infile.is_open()) {fSigCut=3; return;}
214   while(!infile.eof())
215     {
216     infile>>nSigCut;
217   }
218   infile.close();
219   fSigCut=nSigCut; 
220 }//SetSigCutFromFile()    
221 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++  
222 void AliHMPIDCalib::InitHisto(Int_t q,Int_t histocnt,Char_t* name)
223 {
224   //
225   //Init the pad histos. For one DDL we have 11520 pads. ONLY if ENABLED!
226   //Arguments: q-charge, the absolute number of the histogram (AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPcx+1)*(AliHMPIDParam::kMaxPcy+1) and the name of the histogram (unique) 
227   //Returns: none
228   //
229  if(fWritePads==kFALSE) return;
230  fFile->cd();
231  Double_t lowbin,highbin=0;
232  lowbin=q-40.5; highbin=q+40.5;  
233  
234  if(fIsPad[histocnt]==kTRUE) return;
235  
236  if(fLargeHisto==kFALSE) fPadAdc[histocnt]=new TH1I(name,name,81,lowbin,highbin);
237  if(fLargeHisto==kTRUE) fPadAdc[histocnt]=new TH1I(name,name,4093,-0.5,4092.5);
238  fPadAdc[histocnt]->Sumw2();
239  fIsPad[histocnt]=kTRUE;
240  
241 }//InitHisto()
242 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
243 void AliHMPIDCalib::FillHisto(Int_t histocnt,Int_t q)
244 {
245   //
246   //Fill the ADC histograms for each pad
247   //Arguments:  q-charge, the absolute number of the histogram (AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPcx+1)*(AliHMPIDParam::kMaxPcy+1)
248   //Returns: none
249   //
250   if(fIsPad[histocnt]==kFALSE) return;
251   fFile->cd();
252   fPadAdc[histocnt]->Fill(q);
253  
254 }//InitHisto()
255 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
256 void AliHMPIDCalib::InitFile(Int_t ldcId)
257 {
258   //
259   //Initialize the ADC histo output file (one per LDC)
260   //Arguments: LDC Id
261   //Returns: none
262   //
263   if(fWritePads==kFALSE) return;
264   fFile=new TFile(Form("HmpidPadsOnLdc%2d.root",ldcId),"RECREATE");
265 }//InitFile()
266 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
267 void AliHMPIDCalib::CloseFile(Int_t /*ldcId*/)
268 {
269   //
270   //Close the ADC histo output file (one per LDC)
271   //Arguments: LDC Id
272   //Returns: none
273   //
274   fFile->cd();
275   Int_t nPads = (AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPcx+1)*(AliHMPIDParam::kMaxPcy+1);
276   for(Int_t np=0;np<nPads;np++) {if(fIsPad[np]==kTRUE) fPadAdc[np]->Write();} 
277   fFile->Close();
278 }//CloseFile()
279 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
280 void AliHMPIDCalib::FillPedestal(Int_t abspad,Int_t q)
281 {
282   //
283   //Called from the HMPIDda and fills the pedestal values
284   //Arguments: absulote pad number as from AliHMPIDParam and q-charge
285   //Returns: none
286   //
287   if(q<0) AliFatal("Negative charge is read!!!!!!");
288   
289   Int_t nDDL=0, row=0, dil=0, adr=0;
290   //The decoding (abs. pad -> ddl,dil,...) is the same as in AliHMPIDDigit::Raw
291   Int_t y2a[6]={5,3,1,0,2,4};
292
293        nDDL=  2*AliHMPIDParam::A2C(abspad)+AliHMPIDParam::A2P(abspad)%2;              //DDL# 0..13
294  // Int_t tmp=  1+AliHMPIDParam::A2P(abspad)/2*8+AliHMPIDParam::A2Y(abspad)/6;        //temp variable
295   Int_t tmp=   AliHMPIDParam::A2P(abspad)/2*8+AliHMPIDParam::A2Y(abspad)/6;           //temp variable
296 //        row=   (AliHMPIDParam::A2P(abspad)%2)? 25-tmp:tmp;                          //row r=1..24
297         row=   (AliHMPIDParam::A2P(abspad)%2)? tmp:24-tmp;                            //row r=1..24
298         dil=  1+AliHMPIDParam::A2X(abspad)/8;                                         //DILOGIC 
299         adr=y2a[AliHMPIDParam::A2Y(abspad)%6]+6*(AliHMPIDParam::A2X(abspad)%8);       //ADDRESS 0..47 
300   //........... decoding done      
301
302      if(q>0) { 
303         fsq[nDDL][row][dil][adr]+=q;
304       fsq2[nDDL][row][dil][adr]+=q*q;
305       fnpc[nDDL][row][dil][adr]++;                                                     //Count how many times the pad is good (can be different from the good DDL  count)
306                        faddl[nDDL]=kTRUE; 
307                      }
308       else
309       {
310         fpedQ0[nDDL][row][dil][adr]++;                                                 //Count how many times a pad charge is zero
311       }
312       
313      Int_t histocnt=0;   histocnt=(nDDL)*11520+(row-1)*480+(dil-1)*48+adr;             //Histo counter for a single DDL  
314      
315      if(fWritePads==kTRUE)                                                             //works but make it nicer later....
316      { 
317        if( fLargeHisto==kTRUE && nDDL==fSelectDDL) {              
318          InitHisto(q,histocnt,Form("hDDL_%d_Row_%d_Dil_%d_Pad_%d",nDDL,row,dil,adr));  //for large histos use hardware naming
319          FillHisto(histocnt,q);
320         }
321         if(fLargeHisto==kFALSE)
322         {
323          InitHisto(q,histocnt,Form("hPad_Ch_%d_Pc_%d_Px_%d_Py_%d",AliHMPIDParam::A2C(abspad),AliHMPIDParam::A2P(abspad),AliHMPIDParam::A2X(abspad),AliHMPIDParam::A2Y(abspad))); 
324          FillHisto(histocnt,q);  
325         }
326       }//fWritePads
327             
328 }//FillPedestal()
329 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
330 void AliHMPIDCalib::FillErrors(Int_t nDDL,Int_t eType, Int_t nErr)
331 {
332   //
333   //Fill decoding errors from AliHMPIDRawStream
334   //Arguments: nDDL-DDL number, eType- error type as in AliHMPIDRawStream.h and the # of occurence for eType
335   //Retutns: none
336   //
337     if(nErr<=0) return;
338     if(eType < 0 || eType> AliHMPIDRawStream::kSumErr ) return;
339     fErr[nDDL][eType]=fErr[nDDL][eType]+nErr;
340     
341   
342 }//FillErrors()
343 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
344 void AliHMPIDCalib::FillDDLCnt(Int_t iddl,Int_t inDDL, Int_t outDDL)
345 {
346   //
347   //Fill decoding DDL check from RawStream
348   //Arguments: iddl - DDL under setting, inDDL- How many times the DDL is present in the raw stream, outDDL - How many time sthe DDL is succesfylly decoded
349   //Retutns: none
350   //
351  
352   if(inDDL==-1) return;
353   if(fnDDLInStream[iddl]==-1) {fnDDLInStream[iddl]=0; fnDDLOutStream[iddl]=0;}
354   fnDDLInStream[iddl]+=inDDL;
355   fnDDLOutStream[iddl]+=outDDL;
356  
357   
358 }//FillDDLCnt()
359 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
360 Bool_t AliHMPIDCalib::WriteErrors(Int_t nDDL, Char_t* name, Int_t nEv)
361 {
362   //
363   //Write decoding errors to a txt file
364   //Arguments: nDDL-DDL number, name of the error file and number of the read events
365   //Retutns: kTRUE/kFALSE
366   //
367   
368   if(faddl[nDDL]==kFALSE) return kFALSE;                                                                 //if ddl is missing no error file is created
369   ofstream outerr;  outerr.open(name);                                                                   //open error file
370   outerr << Form("%8s %2d\n","RunNumber",(Int_t)fRunNum);                                                //read run number
371   outerr << Form("%8s %2d\n","LdcId" ,          fLdcId);                                                 //read LDC Id
372   outerr << Form("%8s %2d\n","TimeStamp",       fTimeStamp);                                             //read time stamp
373   outerr << Form("%8s %2d\n","TotNumEvt",       nEv);                                                    //read number of total events processed
374   outerr << Form("%8s %2d\n","TotDDLEvt",       fnDDLInStream[nDDL]);                                    //read number of bad events for DDL # nDDL processed
375   outerr << Form("%8s %2d\n","NumBadEvt",       fnDDLInStream[nDDL]-fnDDLOutStream[nDDL]);               //read number of bad events for DDL # nDDL processed
376   outerr << Form("%8s %2.2f\n","NBadE(%)",      (fnDDLInStream[nDDL]-fnDDLOutStream[nDDL])*100.0/nEv);   //read number of bad events (in %) for DDL # nDDL processed
377   
378   for(Int_t  ierr=0; ierr <AliHMPIDRawStream::kSumErr; ierr++) outerr << Form("%2d\t",fErr[nDDL][ierr]); //write errors
379                                                                outerr << Form("\n");                     //last break
380   outerr.close();                                                                                        //write error file
381   return kTRUE;
382     
383 }//FillErrors()
384 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
385 Bool_t AliHMPIDCalib::CalcPedestal(Int_t nDDL, Char_t* name, Int_t nEv)    
386 {
387   //
388   //Calculate pedestal for each pad  
389   //Arguments: nDDL-DDL number, name of the pedestal file and number of the read events
390   //Retutns: kTRUE/kFALSE
391   //
392   
393    
394   if(faddl[nDDL]==kFALSE) return kFALSE;                   //if ddl is missing no ped file is created (and also for LDC selection). Check with Paolo what he checks for?!  
395   Double_t mean=0,sigma=0;
396   Double_t qs2m=0,qsm2=0;
397   ofstream out;                                           //to write the pedestal text files
398   Int_t inhard;
399   Int_t nEvPerPad=0;
400   out.open(name);
401   out << Form("%8s %2d\n","RunNumber",(Int_t)fRunNum);                                                //read run number
402   out << Form("%8s %2d\n","LdcId" ,         fLdcId);                                                  //read LDC Id
403   out << Form("%8s %2d\n","TimeStamp",      fTimeStamp);                                              //read time stamp
404   out << Form("%8s %2d\n","TotNumEvt",      nEv);                                                     //read number of total events processed
405   out << Form("%8s %2d\n","TotDDLEvt",      fnDDLInStream[nDDL]);                                     //read number of bad events for DDL # nDDL processed
406   out << Form("%8s %2d\n","NumBadEvt",      fnDDLInStream[nDDL]-fnDDLOutStream[nDDL]);                //read number of bad events for DDL # nDDL processed
407   out << Form("%8s %2f\n","NBadE(%)",       (fnDDLInStream[nDDL]-fnDDLOutStream[nDDL])*100.0/nEv);    //read number of bad events (in %) for DDL # nDDL processed
408   out << Form("%8s %2.2d\n","#SigCut",      fSigCut);                                                 //# of sigma cuts
409       
410   for(Int_t row = 1; row <= AliHMPIDRawStream::kNRows; row++){
411     for(Int_t dil = 1; dil <= AliHMPIDRawStream::kNDILOGICAdd; dil++){
412       for(Int_t pad = 0; pad < AliHMPIDRawStream::kNPadAdd; pad++){
413         
414         mean  = 50;sigma = 100;
415         
416         nEvPerPad=fnpc[nDDL][row][dil][pad];
417         
418         
419         if(nEvPerPad < 1 ) {                    //if the pad is bad then we assign 100  for the sigma and 50 for the mean
420           mean  = 4000;
421           sigma = 1000;
422         }
423         else{        
424          mean = fsq[nDDL][row][dil][pad]*1.0/nEvPerPad;
425          qs2m = fsq2[nDDL][row][dil][pad]*1.0/nEvPerPad;
426          qsm2 = TMath::Power(fsq[nDDL][row][dil][pad]*1.0/nEvPerPad,2); 
427         sigma = TMath::Sqrt(TMath::Abs(qs2m-qsm2));
428         }
429                 
430         inhard=((Int_t(mean))<<9)+Int_t(mean+3*sigma);
431         out << Form("%2i %2i %2i %5.2f %5.2f %4.4x \n",row,dil,pad,mean,sigma,inhard);
432                  
433        //if(sigma > 3.0) Printf("WARNING SIGMA DDL: %2d row: %2d dil: %2d pad: %2d mean: %3.2f sigma: %2.2f nEvPerPad: %02d fnDDLOutStream: %02d fpedQ0: %02d",nDDL,row,dil,pad,mean,sigma,nEvPerPad,fnDDLOutStream[nDDL],fpedQ0[nDDL][row][dil][pad]);
434         
435        
436         }//adr
437       }//dil
438     }//row
439     out.close();                                          //write pedestal file
440   return kTRUE;
441 }//CaclPedestal()
442 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
443