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