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