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