]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDCalib.cxx
plots (DDL data size and DDL occupancy) scaling corrected
[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 fnDDLInStream(0x0),
48 fnDDLOutStream(0x0),
49 fLargeHisto(kFALSE),
50 fSelectDDL(0),
51 fDeadMap(0x0),
52 fPedMeanMap(0x0),
53 fPedSigMap(0x0),
54 f1DPedMean(0x0),
55 f1DPedSigma(0x0),
56 fNumMaskedPads(0),
57 fNumDeadPads(0)
58 {
59   //
60   //constructor
61   //
62   
63   faddl = new Bool_t[AliHMPIDRawStream::kNDDL];
64   Int_t nPads =  (AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPcx+1)*(AliHMPIDParam::kMaxPcy+1);
65  
66   fpedQ0 = new Int_t***[AliHMPIDRawStream::kNDDL+1];
67   fsq2   = new Float_t ***[AliHMPIDRawStream::kNDDL+1];
68   fsq    = new Float_t ***[AliHMPIDRawStream::kNDDL+1];
69   fnpc   = new Int_t ***[AliHMPIDRawStream::kNDDL+1];
70     fErr = new Int_t*[AliHMPIDRawStream::kNDDL+1];
71    
72   fnDDLInStream  = new Int_t[AliHMPIDRawStream::kNDDL+1];
73   fnDDLOutStream = new Int_t[AliHMPIDRawStream::kNDDL+1];
74
75   
76   for(Int_t iDDL=0;iDDL<AliHMPIDRawStream::kNDDL+1;iDDL++) {
77     
78       fErr[iDDL] = new Int_t[AliHMPIDRawStream::kSumErr+1];
79     fpedQ0[iDDL] = new Int_t**[AliHMPIDRawStream::kNRows+1];
80        fsq[iDDL] = new Float_t**[AliHMPIDRawStream::kNRows+1];
81       fsq2[iDDL] = new Float_t**[AliHMPIDRawStream::kNRows+1];
82       fnpc[iDDL] = new Int_t**[AliHMPIDRawStream::kNRows+1];
83       
84       for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++)  {
85       
86        fpedQ0[iDDL][iRow] = new Int_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
87           fsq[iDDL][iRow] = new Float_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
88          fsq2[iDDL][iRow] = new Float_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
89          fnpc[iDDL][iRow] = new Int_t*[AliHMPIDRawStream::kNDILOGICAdd+1];
90       
91         for(Int_t iDil=1;iDil<AliHMPIDRawStream::kNDILOGICAdd+1;iDil++){
92       
93          fpedQ0[iDDL][iRow][iDil] = new Int_t[AliHMPIDRawStream::kNPadAdd+1];
94            fsq2[iDDL][iRow][iDil] = new Float_t[AliHMPIDRawStream::kNPadAdd+1];
95             fsq[iDDL][iRow][iDil] = new Float_t[AliHMPIDRawStream::kNPadAdd+1];
96            fnpc[iDDL][iRow][iDil] = new Int_t[AliHMPIDRawStream::kNPadAdd+1];
97           }//iDil
98       }//iRow
99    }//iDDL
100     
101    for(Int_t iDDL=0;iDDL<AliHMPIDRawStream::kNDDL+1;iDDL++) {
102         
103      fnDDLInStream[iDDL]=-1;
104      fnDDLOutStream[iDDL]=-1;
105       
106      for(Int_t iErr=0;iErr<AliHMPIDRawStream::kSumErr+1;iErr++)  {fErr[iDDL][iErr]=0;}
107          
108      for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++) {
109         for(Int_t iDil=1;iDil<AliHMPIDRawStream::kNDILOGICAdd+1;iDil++) {
110           for(Int_t iPad=1;iPad<AliHMPIDRawStream::kNPadAdd+1;iPad++) {
111             fpedQ0[iDDL][iRow][iDil][iPad]=0;
112                fsq[iDDL][iRow][iDil][iPad]=0;
113               fsq2[iDDL][iRow][iDil][iPad]=0;
114               fnpc[iDDL][iRow][iDil][iPad]=0;
115         }//iPad
116       }//iDil
117      }//iRow
118    }//iDDL
119     
120   fPadAdc=new TH1I*[nPads];  
121   fIsPad=new Bool_t[nPads];  
122   for(Int_t np=0;np<nPads;np++) {fPadAdc[np]=0x0;   fIsPad[np]=kFALSE;}
123
124   Init();
125 }
126 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
127 AliHMPIDCalib::~AliHMPIDCalib()
128 {
129   //
130   //destructor
131   //
132   if (faddl)     { delete [] faddl;   faddl = 0x0;  } 
133   if (fPadAdc)   { delete [] fPadAdc; fPadAdc=0x0;  }  
134   if (fIsPad)    { delete [] fIsPad;  fIsPad=0x0;   }  
135   if (fFile)     { delete    fFile;   fFile=0x0;    }  
136   if (fDeadMap)  { delete    fDeadMap;fDeadMap=0x0; } 
137   
138   for(Int_t iErr=0;iErr<AliHMPIDRawStream::kSumErr+1;iErr++) { delete [] fErr[iErr];}  delete [] fErr;
139   
140   for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
141    for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++)         
142      for(Int_t iDil=1;iDil<AliHMPIDRawStream::kNDILOGICAdd+1;iDil++)
143       {
144          delete [] fpedQ0[iDDL][iRow][iDil]; //del iPad
145          delete []    fsq[iDDL][iRow][iDil]; //del iPad
146          delete []   fsq2[iDDL][iRow][iDil]; //del iPad
147          delete []   fnpc[iDDL][iRow][iDil]; //del iPad
148        }
149    for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
150      for(Int_t iRow=0;iRow<AliHMPIDRawStream::kNRows+1;iRow++)         
151       {
152         delete [] fpedQ0[iDDL][iRow];  //del iRow
153           delete []  fsq[iDDL][iRow];  //del iRow
154           delete [] fsq2[iDDL][iRow];  //del iRow
155           delete [] fnpc[iDDL][iRow];  //del iRow
156         }
157        
158    for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
159    {   
160        delete [] fpedQ0[iDDL];        //del iRow
161          delete [] fsq2[iDDL];        //del iRow
162          delete []  fsq[iDDL];        //del iRow
163          delete [] fnpc[iDDL];        //del iRow
164      }
165        
166    delete [] fpedQ0;
167    delete [] fsq2;
168    delete [] fsq;
169    delete [] fnpc;
170     
171   fpedQ0=0;    
172     fsq2=0;
173      fsq=0;
174     fnpc=0;
175     
176   fLdcId=0;
177   fTimeStamp=0;
178   fRunNum=0;
179   fSigCut=0;
180   fLargeHisto=kFALSE;
181   fSelectDDL=0;
182   fNumMaskedPads=0;
183   
184   if (fPedMeanMap)   { delete [] fPedMeanMap; fPedMeanMap=0x0;  }  
185   if (fPedSigMap)    { delete [] fPedSigMap;  fPedSigMap=0x0;   }  
186   if (f1DPedMean)    { delete [] f1DPedMean;  f1DPedMean=0x0;   }
187   if (f1DPedSigma)   { delete [] f1DPedSigma; f1DPedSigma=0x0;  }
188   
189 }//dtor
190 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
191 void AliHMPIDCalib::Init()
192 {
193   //
194   //Init the q calc.
195   //Arguments: none
196   //Return: none
197   //
198     
199   fSigCut=3;  //the standard cut
200
201  
202   for(Int_t iDDL=0; iDDL< AliHMPIDRawStream::kNDDL; iDDL++) 
203       {
204          for(Int_t ierr=0; ierr <AliHMPIDRawStream::kSumErr ; ierr++) {
205             fErr[iDDL][ierr]=0;
206             }
207         
208         faddl[iDDL]=kFALSE;
209       }//DDL
210       
211      Int_t      nbins[4]={14, 24, 10, 48};
212      Double_t  binmin[4]={ 0,  1,  1,  0};
213      Double_t  binmax[4]={14, 25, 11, 48};
214      fDeadMap = new THnSparseD("fDeadMap","Dead Channel Map",4,nbins,binmin,binmax); 
215      
216    fPedMeanMap = new TH2F*[AliHMPIDParam::kMaxCh+1];
217    fPedSigMap  = new TH2F*[AliHMPIDParam::kMaxCh+1];
218    f1DPedMean  = new TH1F*[(AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPc+1)];
219    f1DPedSigma = new TH1F*[(AliHMPIDParam::kMaxCh+1)*(AliHMPIDParam::kMaxPc+1)];
220    
221 //   Int_t mapmin,mapmax;
222    for(Int_t iCh=0;iCh<AliHMPIDParam::kMaxCh+1;iCh++)  
223    {
224     fPedMeanMap[iCh]=new TH2F(Form("fPedMeanMap%d",iCh),Form("fPedMeanMap%d;pad x;pad y;Mean pedestal (ADC)",iCh),160, 0,160,144,0,144);fPedMeanMap[iCh]->SetStats(kFALSE);
225     fPedSigMap[iCh] =new TH2F(Form("fPedSigMap%d",iCh), Form("fPedSigMap%d;pad x;pad y;Sigma pedestal (ADC)",iCh), 160, 0,160,144,0,144);fPedSigMap[iCh]->SetStats(kFALSE);
226    }
227    for(Int_t iCh=0;iCh<=AliHMPIDParam::kMaxCh;iCh++)
228    {
229     for(Int_t iFee=0;iFee<6;iFee++)
230      {
231       f1DPedMean[6*iCh+iFee]  = new TH1F(Form("f1DPedMean_Ch%d_FEE_%d" , iCh,iFee),Form("Mean Pedestals, RICH %d, FEE %d;Channels;Mean pedestal (ADC)" ,iCh,iFee),3840,0,3840);f1DPedMean[6*iCh+iFee]->SetStats(kFALSE);
232       f1DPedSigma[6*iCh+iFee] = new TH1F(Form("f1DPedSigma_Ch%d_FEE_%d" ,iCh,iFee),Form("Sigma Pedestal, RICH %d, FEE %d;Channels;Sigma pedestal (ADC)" ,iCh,iFee),3840,0,3840);f1DPedSigma[6*iCh+iFee]->SetStats(kFALSE);
233      } 
234    }
235      
236 }//Init()
237 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
238 void AliHMPIDCalib::SetRunParams(ULong_t runNum,Int_t timeStamp, Int_t ldcId)
239 {
240   //  
241   //Set run parameters for the Pedestal and Error Files
242   //Arguments: run number, time stamp and LDC Id
243   //Returns: none
244   //
245   fRunNum=(Int_t)runNum;
246   fTimeStamp=timeStamp;
247   fLdcId=ldcId;
248 }//SetRunParams()
249 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
250 void AliHMPIDCalib::SetSigCutFromFile(TString hmpInFile)
251 {
252   //
253   //Set Sigma Cut from the file on the LDC, if the input file is not present default value is set!
254   //Arguments: the name of the SigmaCut file on the LDC
255   //Returns: none
256   //
257   Int_t nSigCut=0;
258   ifstream infile(hmpInFile.Data());
259   if(!infile.is_open()) {fSigCut=3; return;}
260   while(!infile.eof())
261     {
262     infile>>nSigCut;
263   }
264   infile.close();
265   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;}
266   Printf("DAQ Sigma Cut from DAQ DB is: %d",nSigCut);
267   fSigCut=nSigCut; 
268 }//SetSigCutFromFile()    
269 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
270 void AliHMPIDCalib::SetDeadChannelMapFromFile(TString hmpInFile)
271 {
272   //
273   //Set Dead Channel Map Cut from the file on the LDC, if the input file is not present default value is set!
274   //Arguments: the name of the Dead Channel Map file on the LDC
275   //Returns: none
276   //
277   Char_t header[256];
278   Int_t ddl=0,row=0,dil=0,pad=0,ch=0,pc=0,chpadx=0,chpady=0,px=0,py=0,isitmasked=0;
279   UInt_t dw=0;
280   Double_t bin[4];
281   ifstream infile(hmpInFile.Data());
282   if(!infile.is_open()) {Printf("HMPID Dead Channel Map file cannot be opened!!!! No mask is applied!!!");return;}
283   infile.getline(header,256);
284   AliHMPIDDigit dig;
285   while(!infile.eof())
286     {
287     infile>>ch>>chpadx>>chpady>>isitmasked;                     //read in masked coordinates; coordinates are in the module coordinate system
288     pc=(chpadx/80)+2*(chpady/48);                               //get PC number
289     px=chpadx-80*(chpadx/80);                                   //get pad X in PC coordinates
290     py=chpady-48*(chpady/48);                                   //get pad Y in PC coordinates --- can we do it better??? -- just with one conversion???? clm
291     if(!dig.Set(ch,pc,px,py,0) && isitmasked) {                 //in the AliHMPIDDigit:Set there is already a check if the coordinates makes sense
292        dig.Raw(dw,ddl,row,dil,pad);
293        bin[0]=ddl; bin[1]=row; bin[2]=dil; bin[3]=pad;
294        fDeadMap->Fill(bin,1);
295       }
296     }
297    infile.close();
298   
299 }//SetDeadChannelMapFromFile()
300 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
301 void AliHMPIDCalib::FillPedestal(Int_t abspad,Int_t q)
302 {
303   //
304   //Called from the HMPIDda and fills the pedestal values
305   //Arguments: absolute pad number as from AliHMPIDParam and q-charge
306   //Returns: none
307   //
308   if(q<0) {
309    AliError("Negative charge is read!!!!!!");
310    return;
311   }
312   UInt_t w32;
313   Int_t nDDL=0, row=0, dil=0, adr=0;
314   //The decoding (abs. pad -> ddl,dil,...) is the same as in AliHMPIDDigit::Raw
315   
316   AliHMPIDDigit dig(abspad,q);
317   dig.Raw(w32,nDDL,row,dil,adr);
318   
319   //........... decoding done      
320
321      if(q>0) { 
322         fsq[nDDL][row][dil][adr]+=q;
323       fsq2[nDDL][row][dil][adr]+=q*q;
324       fnpc[nDDL][row][dil][adr]++;                                                     //Count how many times the pad is good (can be different from the good DDL  count)
325                        faddl[nDDL]=kTRUE; 
326                      }
327       else
328       {
329         fpedQ0[nDDL][row][dil][adr]++;                                                 //Count how many times a pad charge is zero
330       }
331       
332      Int_t histocnt=0;   histocnt=(nDDL)*11520+(row-1)*480+(dil-1)*48+adr;             //Histo counter for a single DDL  
333     
334 }//FillPedestal()
335 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
336 void AliHMPIDCalib::FillErrors(Int_t nDDL,Int_t eType, Int_t nErr)
337 {
338   //
339   //Fill decoding errors from AliHMPIDRawStream
340   //Arguments: nDDL-DDL number, eType- error type as in AliHMPIDRawStream.h and the # of occurence for eType
341   //Retutns: none
342   //
343     if(nErr<=0) return;
344     if(eType < 0 || eType> AliHMPIDRawStream::kSumErr ) return;
345     fErr[nDDL][eType]=fErr[nDDL][eType]+nErr;
346     
347   
348 }//FillErrors()
349 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
350 void AliHMPIDCalib::FillDDLCnt(Int_t iddl,Int_t inDDL, Int_t outDDL)
351 {
352   //
353   //Fill decoding DDL check from RawStream
354   //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
355   //Retutns: none
356   //
357  
358   if(inDDL==-1) return;
359   if(fnDDLInStream[iddl]==-1) {fnDDLInStream[iddl]=0; fnDDLOutStream[iddl]=0;}
360   fnDDLInStream[iddl]+=inDDL;
361   fnDDLOutStream[iddl]+=outDDL;
362  
363   
364 }//FillDDLCnt()
365 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
366 Bool_t AliHMPIDCalib::WriteErrors(Int_t nDDL, Char_t* name, Int_t nEv)
367 {
368   //
369   //Write decoding errors to a txt file
370   //Arguments: nDDL-DDL number, name of the error file and number of the read events
371   //Retutns: kTRUE/kFALSE
372   //  
373   if(faddl[nDDL]==kFALSE) return kFALSE;                                                                 //if ddl is missing no error file is created
374   ofstream outerr;  outerr.open(name);                                                                   //open error file
375   outerr << Form("%8s %2d\n","RunNumber",(Int_t)fRunNum);                                                //read run number
376   outerr << Form("%8s %2d\n","LdcId" ,          fLdcId);                                                 //read LDC Id
377   outerr << Form("%8s %2d\n","TimeStamp",       fTimeStamp);                                             //read time stamp
378   outerr << Form("%8s %2d\n","TotNumEvt",       nEv);                                                    //read number of total events processed
379   outerr << Form("%8s %2d\n","TotDDLEvt",       fnDDLInStream[nDDL]);                                    //read number of bad events for DDL # nDDL processed
380   outerr << Form("%8s %2d\n","NumBadEvt",       fnDDLInStream[nDDL]-fnDDLOutStream[nDDL]);               //read number of bad events for DDL # nDDL processed
381   outerr << Form("%8s %2.2f\n","NBadE(%)",      (fnDDLInStream[nDDL]-fnDDLOutStream[nDDL])*100.0/nEv);   //read number of bad events (in %) for DDL # nDDL processed
382   
383   for(Int_t  ierr=0; ierr <AliHMPIDRawStream::kSumErr; ierr++) outerr << Form("%2d\t",fErr[nDDL][ierr]); //write errors
384                                                                outerr << Form("\n");                     //last break
385   /* write out pads with 0 charge read */
386   for(Int_t row = 1; row <= AliHMPIDRawStream::kNRows; row++){
387     for(Int_t dil = 1; dil <= AliHMPIDRawStream::kNDILOGICAdd; dil++){
388       for(Int_t pad = 0; pad < AliHMPIDRawStream::kNPadAdd; pad++){
389         if(fpedQ0[nDDL][row][dil][pad]>0) outerr<< Form("%2d %2d %2d %3d\n",row,dil,pad,fpedQ0[nDDL][row][dil][pad]);
390       }
391     }
392   } 
393                                                                                                                                                                                        
394                                                                
395   outerr.close();                                                                                        //write error file
396   
397   return kTRUE;
398     
399 }//FillErrors()
400 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
401 Bool_t AliHMPIDCalib::CalcPedestal(Int_t nDDL, Char_t* name, Char_t *name2,Int_t nEv)    
402 {
403   //
404   //Calculate pedestal for each pad  
405   //Arguments: nDDL-DDL number, name of the pedestal file and number of the read events
406   //Retutns: kTRUE/kFALSE
407   //
408   
409   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?!  
410
411   Int_t feeOffset=196657;
412   ofstream feeInput; feeInput.open(Form("%s",name2));      //write thr file for Fe2C
413   
414   Double_t mean=0,sigma=0, threshold=0;
415   Double_t qs2m=0,qsm2=0;
416   ofstream out;                                            //to write the pedestal text files
417   Int_t inhard;
418   Int_t nEvPerPad=0;
419   Int_t pedbin=0;
420   
421   Int_t abspad,ch,pc,pcx,pcy,chX,chY,fee;
422   Int_t binSp[4]={0};
423   out.open(name);
424   out << Form("%8s %2d\n","RunNumber",(Int_t)fRunNum);                                                //read run number
425   out << Form("%8s %2d\n","LdcId" ,         fLdcId);                                                  //read LDC Id
426   out << Form("%8s %2d\n","TimeStamp",      fTimeStamp);                                              //read time stamp
427   out << Form("%8s %2d\n","TotNumEvt",      nEv);                                                     //read number of total events processed
428   out << Form("%8s %2d\n","TotDDLEvt",      fnDDLInStream[nDDL]);                                     //read number of bad events for DDL # nDDL processed
429   out << Form("%8s %2d\n","NumBadEvt",      fnDDLInStream[nDDL]-fnDDLOutStream[nDDL]);                //read number of bad events for DDL # nDDL processed
430   out << Form("%8s %2f\n","NBadE(%)",       (fnDDLInStream[nDDL]-fnDDLOutStream[nDDL])*100.0/nEv);    //read number of bad events (in %) for DDL # nDDL processed
431   out << Form("%8s %d\n","#SigCut",      fSigCut);                                                    //# of sigma cuts
432       
433   for(Int_t row = 1; row <= AliHMPIDRawStream::kNRows; row++){
434     feeInput << Form("0xabcdabcd \n");                                                                    //before each row we write a marker to separate the rows within a DDL                       
435     
436    
437     for(Int_t dil = 1; dil <= AliHMPIDRawStream::kNDILOGICAdd; dil++){
438       for(Int_t pad = 0; pad < AliHMPIDRawStream::kNPadAdd; pad++){
439         mean  = 50;sigma = 100;                                                   //init maen and sigma to a low value
440         nEvPerPad=fnpc[nDDL][row][dil][pad];                                      //check how many times the pad was read out
441         abspad=AliHMPIDRawStream::GetPad(nDDL,row,dil,pad);                       //get the absolute oad coordinate
442         ch=AliHMPIDParam::A2C(abspad);                                            //get chamber number
443         pc=AliHMPIDParam::A2P(abspad);                                            //get PC number
444         pcx=AliHMPIDParam::A2X(abspad);                                           //get pad x in PC
445         pcy=AliHMPIDParam::A2Y(abspad);                                           //get pad y in PC
446         chX = (pc%2)*AliHMPIDParam::kPadPcX+pcx;                                  //get pad x in Ch   
447         chY = (pc/2)*AliHMPIDParam::kPadPcY+pcy;                                  //get pad y in Ch
448         binSp[0]=nDDL+1;binSp[1]=row;binSp[2]=dil;binSp[3]=pad+1;                 //set dead map coordinates for check
449         
450        if(nEvPerPad < 1 ) {                                                      //if the pad is bad then we assign 100  for the sigma and 50 for the mean
451           mean  = AliHMPIDParam::kPadMeanZeroCharge;
452           sigma = AliHMPIDParam::kPadSigmaZeroCharge;
453           fNumDeadPads++;
454         }
455         else if(fDeadMap->GetBinContent(binSp)>0)                                 //check if channel is masked, if yes set maksed values
456         {
457           mean  = AliHMPIDParam::kPadMeanMasked;
458           sigma = AliHMPIDParam::kPadSigmaMasked;
459           fNumMaskedPads++;
460         }
461        else{            
462          mean = fsq[nDDL][row][dil][pad]*1.0/nEvPerPad;
463          qs2m = fsq2[nDDL][row][dil][pad]*1.0/nEvPerPad;
464          qsm2 = TMath::Power(fsq[nDDL][row][dil][pad]*1.0/nEvPerPad,2); 
465         sigma = TMath::Sqrt(TMath::Abs(qs2m-qsm2));
466         }
467         
468         //The electronics takes the 32bit int as: first 9 bits for the pedestal and the second 9 bits for threshold
469         threshold = mean+fSigCut*sigma;                                                                    
470         if(mean > 511.0 || threshold > 511.0) {mean = AliHMPIDParam::kPadMeanMasked; threshold = AliHMPIDParam::kPadMeanMasked + 5.0 * AliHMPIDParam::kPadSigmaMasked; }
471         //inhard=((Int_t(mean+fSigCut*sigma))<<9)+Int_t(mean);                                                 //right calculation, xchecked with Paolo 8/4/2008
472         inhard=((Int_t(threshold))<<9)+Int_t(mean);                                                            //right calculation, xchecked with Paolo 8/4/2008
473         
474         out << Form("%2i %2i %2i %5.3f %5.3f %4.4x \n",row,dil,pad,mean,sigma,inhard);
475         feeInput << Form("0x%4.4x\n",inhard);
476         
477         // fill histograms to be exported to AMORE    
478         fPedMeanMap[ch]->SetTitle(Form("PedMeanMap%d RunNum: %d",ch,fRunNum));
479         fPedSigMap[ch]->SetTitle(Form("PedSigmaMap%d RunNum: %d",ch,fRunNum));
480         fPedMeanMap[ch]->Fill(chX,chY,mean);         
481         fPedSigMap[ch]->Fill(chX,chY,sigma);         
482         if(nDDL%2==0) pedbin = (24-row)*2*480+(10-dil)*48+pad;
483         if(nDDL%2!=0) pedbin = (row*2-1)*480+(10-dil)*48+pad;
484         pedbin = pedbin - 3840*(pedbin/3840);
485         fee=AliHMPIDRawStream::GetFee(nDDL,row);
486         f1DPedMean[6*(nDDL/2)+fee]->SetTitle(Form("PedMean_Ch%d_FEE_%d RunNum: %d",ch,fee,fRunNum));
487         f1DPedSigma[6*(nDDL/2)+fee]->SetTitle(Form("PedSigma_Ch%d_FEE_%d RunNum: %d",ch,fee,fRunNum));
488         f1DPedMean[6*(nDDL/2)+fee]->Fill(pedbin,mean);
489         f1DPedSigma[6*(nDDL/2)+fee]->Fill(pedbin,sigma);
490
491         }//adr==pad
492         //we have to write up to 64 not 48 in the DILOGIC since they are daisy chained!
493         //offset and format is defined for the Fe2C code
494         for(Int_t idd=0;idd<16;idd++) feeInput << Form("0x%4.4x\n",idd+feeOffset);                 
495       }//dil
496       
497       
498     }//row
499     out.close();                                          //write pedestal file
500     feeInput.close();
501  
502   return kTRUE;
503 }//CaclPedestal()
504 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++