]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDPreprocessor.cxx
Coverity fix
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDPreprocessor.cxx
1 #include "AliHMPIDPreprocessor.h" //header no includes
2 #include "AliHMPIDDigit.h"        //ProcPed()
3 #include "AliHMPIDRawStream.h"    //ProcPed()
4 #include <Riostream.h>            //ProcPed()  
5 #include <AliLog.h>               //all
6 #include <AliCDBMetaData.h>       //ProcPed(), ProcDcs()
7 #include <AliDCSValue.h>          //ProcDcs()
8 #include <TObjString.h>           //ProcDcs(), ProcPed()
9 #include <TTimeStamp.h>           //Initialize()
10 #include <TF1.h>                  //Process()
11 #include <TF2.h>                  //Process()
12 #include <TString.h>
13 #include <TGraph.h>               //Process()
14 #include <TMatrix.h>              //ProcPed()
15 #include <TList.h>                //ProcPed()
16 #include <TSystem.h>              //ProcPed()
17 //.
18 // HMPID Preprocessor base class
19 //.
20 //.
21 //.
22 ClassImp(AliHMPIDPreprocessor)
23
24 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
25 void AliHMPIDPreprocessor::Initialize(Int_t run, UInt_t startTime,UInt_t endTime)
26 {
27 // Initialize the parameter coming from AliPreprocessor
28 //  run -> run number
29 // startTime -> starting time 
30 // endTime   -> ending time
31   AliPreprocessor::Initialize(run, startTime, endTime);
32   
33   AliInfo(Form("HMPID started for Run %d \n\tStartTime %s \n\t  EndTime %s", run,TTimeStamp(startTime).AsString(),TTimeStamp(endTime).AsString()));
34
35 }
36 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
37 UInt_t AliHMPIDPreprocessor::Process(TMap* pMap)
38 {
39 // Process all information from DCS and DAQ
40 // Arguments: pMap- map of DCS aliases
41 // Returns: 0 on success or 1 on error (opposite to Store!)
42
43   TString runType = GetRunType();
44   Log(Form(" AliHMPIDPreprocessor: RunType is %s",runType.Data()));
45   Bool_t statusNoise=kFALSE, statusDcs=kFALSE;
46 // start to check event type and procedures
47   
48   Log("HMPID - Process in Preprocessor started");
49   if(! pMap) {
50     Log("HMPID - ERROR - Not map of DCS aliases for HMPID - ");             return kTRUE;   // error in the DCS mapped aliases
51   }   
52   if (runType == "CALIBRATION"){
53     if (!ProcPed()){
54         Log("HMPID - ERROR - Pedestal processing failed!!");                return kTRUE;   // error in pedestal processing
55     } else {
56         Log("HMPID - Pedestal processing successful!!");                    return kFALSE;  // ok for pedestals
57     }
58   }//CALIBRATION
59   else if ( runType=="STANDALONE" || runType=="PHYSICS"){   
60    statusDcs=ProcDcs(pMap);
61    statusNoise=ProcNoiseMap();
62    if(!statusDcs || !statusNoise) { Log(Form("HMPID - ERROR - Noise Map(%d) and/or DCS(%d) processing failed!! (0=OK, 1=FAILED)",statusNoise,statusDcs)); return kTRUE; }  // error in Noise Map or DCS processing
63    else                           { Log("HMPID - Noise Map and DCS processing successful!!");  return kFALSE;}  // ok
64   }//STANDALONE or PHYSICS run
65    else {
66     Log("HMPID - Nothing to do with preprocessor for HMPID, bye!");         return kFALSE;  // ok - nothing done
67   }
68 }//Process()
69 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
70 Bool_t AliHMPIDPreprocessor::ProcNoiseMap()
71 {
72   //
73   // Goal: Process the Noise Map created by the HMP Physics DA to mask 
74   // noisy channels from reconstruction and follow changes in accepatnce
75   // eg. DDL turn on/off after PEDESTAL run and between PHYSICS runs.
76   // Returns kFALSE on success
77  
78   Bool_t stProcNoise=kFALSE;
79   TFile  *fNoiseFile;
80   TH2F   *hNoiseMap = 0x0;
81   
82   TList *pNoiseSource=GetFileSources(kDAQ,"HmpPhysicsDaNoiseMap.root"); //get list of DAQ source names containing id "HmpPhysicsDaNoiseMap" --> defined in HMPIDphysda.cxx
83   if(!pNoiseSource) {Log(Form("ERROR: Retrieval of sources for noise map: HmpPhysicsDaNoiseMap.root is failed!")); return stProcNoise;}
84   if(!(TObjString*)pNoiseSource->At(0)) {Log(Form("ERROR: empty list received from DAQ Source!")); return stProcNoise;}
85     
86   TString noiseFile = GetFile(kDAQ,Form("HmpPhysicsDaNoiseMap.root"),((TObjString*)pNoiseSource->At(0))->GetName());
87   if(noiseFile.Length()==0) {Log(Form("ERROR retrieving noise map file: HmpPhysicsDaNoiseMap.root")); return stProcNoise;}
88   
89   fNoiseFile = TFile::Open(noiseFile.Data(),"read");
90   if(!fNoiseFile) {Log(Form("ERROR cannot open NoiseFile: %s!",noiseFile.Data())); return stProcNoise;}
91   hNoiseMap = (TH2F*) fNoiseFile->Get("hHmpNoiseMaps");
92   
93   AliCDBMetaData metaDataHisto;
94   metaDataHisto.SetBeamPeriod(0);
95   metaDataHisto.SetResponsible("AliHMPIDPreprocessor"); 
96   metaDataHisto.SetComment("AliHMPIDPreprocessor stores the Noise Map object as Reference Data.");
97   AliInfo("Storing Reference Data");
98   stProcNoise = Store("Calib","NoiseMap",hNoiseMap,&metaDataHisto,0,kTRUE);
99   if(!stProcNoise) {
100     Log("HMPID - failure to store Noise Map data results in OCDB");   
101   }
102   return stProcNoise;
103 }//ProcNoiseMap
104 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
105 Bool_t AliHMPIDPreprocessor::ProcDcs(TMap* pMap)
106 {
107 // Process: 1 (old). inlet and outlet C6F14 temperature, stores TObjArray of 21 TF1, where TF1 is Nmean=f(t), one per radiator
108 // Process: 1. inlet and outlet C6F14 temperature, stores TObjArray of 42 TF1, where TF1 are Tin and Tout per radiator
109 //             + one function for the mean energy photon (in total 43).
110 //          2. CH4 pressure and HV                 stores TObjArray of 7 TF1 where TF1 is thr=f(t), one per chamber
111 // Arguments: pDcsMap - map of structure "alias name" - TObjArray of AliDCSValue
112 // Assume that: HV is the same during the run for a given chamber, different chambers might have different HV
113 //              P=f(t), different for different chambers
114 // Returns: kTRUE on success  
115
116   Bool_t stDcsStore=kFALSE;
117
118 // Qthr=f(HV,P) [V,mBar]  logA0=k*HV+b is taken from p. 64 TDR plot 2.59 for PC32 
119 //                           A0=f(P) is taken from DiMauro mail
120 // Qthr is estimated as 3*A0
121
122   TF2 thr("RthrCH4"  ,"3*10^(3.01e-3*x-4.72)+170745848*exp(-y*0.0162012)"             ,2000,3000,900,1200); 
123   
124   TObjArray arNmean(43);       arNmean.SetOwner(kTRUE);     //42 Tin and Tout one per radiator + 1 for ePhotMean
125   TObjArray arQthre(42);       arQthre.SetOwner(kTRUE);     //42 Qthre=f(time) one per sector
126   
127   AliDCSValue *pVal; Int_t cnt=0;
128   
129   Double_t xP,yP;
130
131 //  TF1 **pTin  = new TF1*[21];
132 //  TF1 **pTout = new TF1*[21];
133   TF1  *pTin[21];
134   TF1 *pTout[21];
135
136 // evaluate Environment Pressure
137   
138   TObjArray *pPenv=(TObjArray*)pMap->GetValue("HMP_DET/HMP_ENV/HMP_ENV_PENV.actual.value");
139   if(!pPenv) {
140     AliWarning(" No Data Points from HMP_ENV_PENV.actual.value!");
141     return kFALSE;
142   } else {
143     Log(Form(" Environment Pressure data              ---> %3i entries",pPenv->GetEntries()));
144     if(pPenv->GetEntries()) {
145       TIter nextPenv(pPenv);
146       TGraph *pGrPenv=new TGraph; cnt=0;
147       while((pVal=(AliDCSValue*)nextPenv())) pGrPenv->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());        //P env
148       if( cnt==1) {
149         pGrPenv->GetPoint(0,xP,yP);
150         new TF1("Penv",Form("%f",yP),fStartTime,fEndTime);
151       } else {
152         pGrPenv->Fit(new TF1("Penv","1000+x*[0]",fStartTime,fEndTime),"Q");
153       }
154       delete pGrPenv;
155     } else {AliWarning(" No Data Points from HMP_ENV_PENV.actual.value!");return kFALSE;}
156   }  
157 // evaluate Pressure
158   
159   for(Int_t iCh=0;iCh<7;iCh++){                   
160     TObjArray *pP =(TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_GAS/HMP_MP%i_GAS_PMWPC.actual.value",iCh,iCh,iCh));
161     if(!pP) {
162       AliWarning(Form(" No Data Points from HMP_MP%1i_GAS_PMWPC.actual.value!",iCh));
163       return kFALSE;
164     } else {
165         Log(Form(" Pressure for module %i data             ---> %3i entries",iCh,pP->GetEntries()));
166       if(pP->GetEntries()) {
167         TIter nextP(pP);    
168         TGraph *pGrP=new TGraph; cnt=0; 
169         while((pVal=(AliDCSValue*)nextP())) pGrP->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());            //P
170         if( cnt==1) {
171           pGrP->GetPoint(0,xP,yP);
172           new TF1(Form("P%i",iCh),Form("%f",yP),fStartTime,fEndTime);
173         } else {
174           pGrP->Fit(new TF1(Form("P%i",iCh),"[0] + x*[1]",fStartTime,fEndTime),"Q");
175         }
176         delete pGrP;
177       } else {AliWarning(" No Data Points from HMP_MP0-6_GAS_PMWPC.actual.value!");return kFALSE;}
178     }
179     
180 // evaluate High Voltage
181     
182    for(Int_t iSec=0;iSec<6;iSec++){
183       TObjArray *pHV=(TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_PW/HMP_MP%i_SEC%i/HMP_MP%i_SEC%i_HV.actual.vMon",iCh,iCh,iCh,iSec,iCh,iSec));
184       if(!pHV) {
185         AliWarning(Form(" No Data Points from HMP_MP%1i_SEC%1i_HV.actual.vMon!",iCh,iSec));
186         return kFALSE;
187       } else {
188         Log(Form(" HV for module %i and secto %i data       ---> %3i entries",iCh,iSec,pHV->GetEntries()));
189         if(pHV->GetEntries()) {
190           TIter nextHV(pHV);
191           TGraph *pGrHV=new TGraph; cnt=0;
192           while((pVal=(AliDCSValue*)nextHV())) pGrHV->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat());            //HV
193           if( cnt==1) {
194             pGrHV->GetPoint(0,xP,yP);
195             new TF1(Form("HV%i_%i",iCh,iSec),Form("%f",yP),fStartTime,fEndTime);               
196           } else {
197             pGrHV->Fit(new TF1(Form("HV%i_%i",iCh,iSec),"[0]+x*[1]",fStartTime,fEndTime),"Q");               
198           }
199           delete pGrHV;
200         } else {AliWarning(" No Data Points from HMP_MP0-6_SEC0-5_HV.actual.vMon!");return kFALSE;}
201       }   
202 // evaluate Qthre
203      
204       arQthre.AddAt(new TF1(Form("HMP_QthreC%iS%i",iCh,iSec),
205           Form("3*10^(3.01e-3*HV%i_%i - 4.72)+170745848*exp(-(P%i+Penv)*0.0162012)",iCh,iSec,iCh),fStartTime,fEndTime),6*iCh+iSec);
206
207      //arQthre.AddAt(new TF1(Form("HMP_QthreC%iS%i",iCh,iSec),"100",fStartTime,fEndTime),6*iCh+iSec);  
208     }
209 // evaluate Temperatures: in and out of the radiators    
210     // T in
211     for(Int_t iRad=0;iRad<3;iRad++){
212       
213       pTin[3*iCh+iRad]  = new TF1(Form("Tin%i%i" ,iCh,iRad),"[0]+[1]*x",fStartTime,fEndTime);
214       pTout[3*iCh+iRad] = new TF1(Form("Tout%i%i",iCh,iRad),"[0]+[1]*x",fStartTime,fEndTime);          
215
216       //pTin[3*iCh+iRad]  = new TF1(Form("Tin%i%i" ,iCh,iRad),"21",fStartTime,fEndTime);
217       //pTout[3*iCh+iRad] = new TF1(Form("Tout%i%i",iCh,iRad),"22",fStartTime,fEndTime);          
218                   
219     TObjArray *pT1=(TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_LIQ_LOOP.actual.sensors.Rad%iIn_Temp",iCh,iCh,iRad));
220       if(!pT1) {
221         AliWarning(Form(" No Data Points from HMP_MP%1i_LIQ_LOOP.actual.sensors.Rad%1iIn_Temp!",iCh,iRad));
222         return kFALSE;
223       } else {
224         Log(Form(" Temperatures for module %i inside data  ---> %3i entries",iCh,pT1->GetEntries()));
225         if(pT1->GetEntries()) {
226           TIter nextT1(pT1);//Tin
227           TGraph *pGrT1=new TGraph; cnt=0;  while((pVal=(AliDCSValue*)nextT1())) pGrT1->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat()); //T inlet
228           if(cnt==1) { 
229             pGrT1->GetPoint(0,xP,yP);
230             pTin[3*iCh+iRad]->SetParameter(0,yP);
231             pTin[3*iCh+iRad]->SetParameter(1,0);
232           } else {
233             pGrT1->Fit(pTin[3*iCh+iRad],"Q");
234           }
235           delete pGrT1;
236         } else {AliWarning(" No Data Points from HMP_MP0-6_LIQ_LOOP.actual.sensors.Rad0-2In_Temp!");return kFALSE;}
237       }
238     // T out
239       TObjArray *pT2=(TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_MP%i/HMP_MP%i_LIQ_LOOP.actual.sensors.Rad%iOut_Temp",iCh,iCh,iRad)); 
240       if(!pT2) {
241         AliWarning(Form(" No Data Points from HMP_MP%1i_LIQ_LOOP.actual.sensors.Rad%1iOut_Temp!",iCh,iRad));
242         return kFALSE;
243       } else {
244         Log(Form(" Temperatures for module %i outside data ---> %3i entries",iCh,pT2->GetEntries()));
245         if(pT2->GetEntries()) {
246           TIter nextT2(pT2);//Tout      
247           TGraph *pGrT2=new TGraph; cnt=0;  while((pVal=(AliDCSValue*)nextT2())) pGrT2->SetPoint(cnt++,pVal->GetTimeStamp(),pVal->GetFloat()); //T outlet 
248           if(cnt==1) { 
249             pGrT2->GetPoint(0,xP,yP);
250             pTout[3*iCh+iRad]->SetParameter(0,yP);
251             pTout[3*iCh+iRad]->SetParameter(1,0);
252           } else {
253             pGrT2->Fit(pTout[3*iCh+iRad],"Q");
254           }
255           delete pGrT2;
256         } else {AliWarning(" No Data Points from HMP_MP0-6_LIQ_LOOP.actual.sensors.Rad0-2Out_Temp!");return kFALSE;}
257       }
258         
259 // evaluate Mean Refractive Index
260       
261       arNmean.AddAt(pTin[3*iCh+iRad] ,6*iCh+2*iRad  ); //Tin =f(t)
262       arNmean.AddAt(pTout[3*iCh+iRad],6*iCh+2*iRad+1); //Tout=f(t)
263       
264     }//radiators loop
265   }//chambers loop
266   
267   Double_t eMean = ProcTrans(pMap);
268   arNmean.AddAt(new TF1("HMP_PhotEmean",Form("%f",eMean),fStartTime,fEndTime),42); //Photon energy mean
269     
270   AliCDBMetaData metaData; 
271   metaData.SetBeamPeriod(0); 
272   metaData.SetResponsible("AliHMPIDPreprocessor"); 
273   metaData.SetComment("HMPID preprocessor fills TObjArrays.");
274
275   stDcsStore =   Store("Calib","Qthre",&arQthre,&metaData,0,kTRUE) &&    // from DCS  0,kTRUE generates the file from Run 0 to Run 99999999
276                  Store("Calib","Nmean",&arNmean,&metaData,0,kTRUE);      // from DCS
277 //  stDcsStore =   Store("Calib","Qthre",&arQthre,&metaData) &&    // from DCS 
278 //                 Store("Calib","Nmean",&arNmean,&metaData);      // from DCS
279   if(!stDcsStore) {
280     Log("HMPID - failure to store DCS data results in OCDB");    
281   }
282   
283 //  arNmean.Delete();
284 //  arQthre.Delete();
285
286 //  for(Int_t i=0;i<21;i++) delete pTin[i]; delete []pTin;
287 //  for(Int_t i=0;i<21;i++) delete pTout[i]; delete []pTout;
288   
289   return stDcsStore;
290 }//Process()
291 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
292 Bool_t AliHMPIDPreprocessor::ProcPed()
293 {
294 // Process pedestal files and create 7 M(padx,pady)=sigma, one for each chamber
295 // Arguments:
296 // Returns: kTRUE on success
297   
298   Bool_t stPedStore=kFALSE;
299   Bool_t stDeadMaskedStore=kFALSE;
300   AliHMPIDDigit dig;
301   AliHMPIDRawStream rs;
302   Int_t nSigCut,r,d,a,hard;  Float_t mean,sigma;
303   Int_t  runNumber,ldcId,timeStamp,nEv,nDdlEv,nBadEv;  Char_t tName[10]; 
304   Float_t nBadEvPer;
305   
306   TObjArray aDaqSig(7);     aDaqSig.SetOwner(kTRUE);     for(Int_t i=0;i<7;i++) aDaqSig.AddAt(new TMatrix(160,144),i);         //TObjArray of 7 TMatrixF, m(padx,pady)=sigma
307   TObjArray aDeadMasked(7); aDeadMasked.SetOwner(kTRUE); for(Int_t i=0;i<7;i++) aDeadMasked.AddAt(new TMatrix(160,144),i);     //TObjArray of 7 TMatrixF, m(padx,pady)=pedestal
308   
309   for(Int_t iddl=0;iddl<AliHMPIDRawStream::kNDDL;iddl++)            //retrieve the files from LDCs independently the DDL<->LDC connection
310   {
311     TList *pLdc=GetFileSources(kDAQ,Form("HmpidPedDdl%02i.txt",iddl)); //get list of LDC names containing id "pedestals"
312     if(!pLdc) {Log(Form("ERROR: Retrieval of sources for pedestals: HmpidPedDdl%02i.txt failed!",iddl));continue;}
313     
314     Log(Form("HMPID - Pedestal files to be read --> %i LDCs for HMPID",pLdc->GetEntries()));
315     for(Int_t i=0;i<pLdc->GetEntries();i++) {//lists of LDCs -- but in general we have 1 LDC for 1 ped file
316     TString fileName = GetFile(kDAQ,Form("HmpidPedDdl%02i.txt",iddl),((TObjString*)pLdc->At(i))->GetName());
317     if(fileName.Length()==0) {Log(Form("ERROR retrieving pedestal file: HmpidPedDdl%02i.txt!",iddl));continue;}
318   
319     //reading pedestal file
320     ifstream infile(fileName.Data()); 
321     
322     if(!infile.is_open()) {Log("No pedestal file found for HMPID,bye!");continue;}
323     TMatrix *pM=(TMatrixF*)aDaqSig.At(iddl/2);
324     TMatrix *pDM=(TMatrixF*)aDeadMasked.At(iddl/2);
325     infile>>tName>>runNumber;Printf("Xcheck: reading run %i",runNumber);
326     infile>>tName>>ldcId;
327     infile>>tName>>timeStamp;
328     infile>>tName>>nEv; 
329     infile>>tName>>nDdlEv;
330     infile>>tName>>nBadEv;
331     infile>>tName>>nBadEvPer;
332     infile>>tName>>nSigCut; pM->SetUniqueID(nSigCut); //n. of pedestal distribution sigmas used to create zero suppresion table
333     while(!infile.eof()){
334       infile>>dec>>r>>d>>a>>mean>>sigma>>hex>>hard;
335       if(rs.GetPad(iddl,r,d,a)>=0){                     //the GetPad returns meaningful abs pad number                                                          
336       dig.SetPad(rs.GetPad(iddl,r,d,a));
337       dig.SetQ((Int_t)mean);
338       (*pM)(dig.PadChX(),dig.PadChY()) = sigma;
339       if( (mean == AliHMPIDParam::kPadMeanZeroCharge && sigma == AliHMPIDParam::kPadSigmaZeroCharge) || 
340           (mean == AliHMPIDParam::kPadMeanMasked     && sigma == AliHMPIDParam::kPadSigmaMasked)     ) 
341         {(*pDM)(dig.PadChX(),dig.PadChY()) = mean;} 
342       }
343     }
344     infile.close();
345     Log(Form("Pedestal file for DDL %i read successfully",iddl));
346   
347  
348     
349   }//LDCs reading entries
350
351  }//DDL 
352
353   AliCDBMetaData metaData; 
354   metaData.SetBeamPeriod(0); 
355   metaData.SetResponsible("AliHMPIDPreprocessor"); 
356   metaData.SetComment("HMPID processor fills TObjArrays.");  
357   stPedStore = Store("Calib","DaqSig",&aDaqSig,&metaData,0,kTRUE);
358   if(!stPedStore) {     Log("HMPID - failure to store PEDESTAL data results in OCDB");      }
359   stDeadMaskedStore = Store("Calib","Masked",&aDeadMasked,&metaData,0,kTRUE);
360   if(!stDeadMaskedStore) {     Log("HMPID - failure to store DEAD & MASKED channel map  in OCDB");      }
361   Bool_t pedRes=stPedStore*stDeadMaskedStore;
362   return pedRes;
363   
364 }//ProcPed()  
365 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
366 Double_t AliHMPIDPreprocessor::ProcTrans(TMap* pMap)
367 {
368   //  Process transparency monitoring data and calculates Emean  
369   
370   Double_t sEnergProb=0, sProb=0;
371
372   Double_t tRefCR5 = 19. ;                                      // mean temperature of CR5 where the system is in place
373
374   Double_t eMean = 0;
375       
376   AliDCSValue *pVal;
377
378   for(Int_t i=0; i<30; i++){
379
380     // evaluate wavelenght 
381     TObjArray *pWaveLenght = (TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.waveLenght",i));
382     if(!pWaveLenght){ 
383         AliWarning(Form("No Data Point values for HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.waveLenght  -----> Default E mean used!!!!!",i));
384         return DefaultEMean(); // to be checked
385       } 
386         
387     pVal=(AliDCSValue*)pWaveLenght->At(0);
388     Double_t lambda = pVal->GetFloat();
389
390     if(lambda<150. || lambda>230.){ 
391         AliWarning(Form("Wrong value for HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.waveLenght  -----> Default E mean used!!!!!",i));
392         return DefaultEMean(); // to be checked
393       } 
394
395     Double_t photEn = 1239.842609/lambda;     // 1239.842609 from nm to eV
396     
397     if(photEn<AliHMPIDParam::EPhotMin() || photEn>AliHMPIDParam::EPhotMax()) continue;
398     
399     // evaluate phototube current for argon reference
400     TObjArray *pArgonRef  = (TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.argonReference",i));
401     if(!pArgonRef){ 
402         AliWarning(Form("No Data Point values for HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.argonReference  -----> Default E mean used!!!!!",i));
403         return DefaultEMean(); // to be checked
404       } 
405
406     pVal=(AliDCSValue*)pArgonRef->At(0);    
407     Double_t aRefArgon = pVal->GetFloat();
408
409     // evaluate phototube current for argon cell
410     TObjArray *pArgonCell = (TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.argonCell",i));
411     if(!pArgonCell){ 
412         AliWarning(Form("No Data Point values for HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.argonCell  -----> Default E mean used!!!!!",i));
413         return DefaultEMean(); // to be checked
414       } 
415       
416     pVal=(AliDCSValue*)pArgonRef->At(0);
417     Double_t aCellArgon = pVal->GetFloat();
418
419     //evaluate phototube current for freon reference
420     TObjArray *pFreonRef  = (TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.c6f14Reference",i));
421     if(!pFreonRef){ 
422         AliWarning(Form("No Data Point values for HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.c6f14Reference  -----> Default E mean used!!!!!",i));
423         return DefaultEMean(); // to be checked
424       } 
425         
426     pVal=(AliDCSValue*)pFreonRef->At(0);
427     Double_t aRefFreon = pVal->GetFloat();
428
429     //evaluate phototube current for freon cell
430     TObjArray *pFreonCell = (TObjArray*)pMap->GetValue(Form("HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.c6f14Cell",i));
431     if(!pFreonCell){
432         AliWarning(Form("No Data Point values for HMP_DET/HMP_INFR/HMP_INFR_TRANPLANT/HMP_INFR_TRANPLANT_MEASURE.mesure%i.c6f14Cell  -----> Default E mean used!!!!!",i));
433         return DefaultEMean(); // to be checked
434       }
435
436     pVal=(AliDCSValue*)pFreonCell->At(0);
437     Double_t aCellFreon = pVal->GetFloat();
438  
439    //evaluate correction factor to calculate trasparency (Ref. NIMA 486 (2002) 590-609)
440     
441     Double_t aN1 = AliHMPIDParam::NIdxRad(photEn,tRefCR5);
442     Double_t aN2 = AliHMPIDParam::NMgF2Idx(photEn);
443     Double_t aN3 = 1;                              // Argon Idx
444
445     Double_t aR1               = ((aN1 - aN2)*(aN1 - aN2))/((aN1 + aN2)*(aN1 + aN2));
446     Double_t aR2               = ((aN2 - aN3)*(aN2 - aN3))/((aN2 + aN3)*(aN2 + aN3));
447     Double_t aT1               = (1 - aR1);
448     Double_t aT2               = (1 - aR2);
449     Double_t aCorrFactor       = (aT1*aT1)/(aT2*aT2);
450
451     // evaluate 15 mm of thickness C6F14 Trans
452     Double_t aTransRad;
453     
454     Double_t aConvFactor = 1.0 - 0.3/1.8;         
455         
456     if(aRefFreon*aRefArgon>0) {
457       aTransRad  = TMath::Power((aCellFreon/aRefFreon)/(aCellArgon/aRefArgon)*aCorrFactor,aConvFactor);
458     } else {
459       return DefaultEMean();
460     }
461
462     // evaluate 0.5 mm of thickness SiO2 Trans
463     Double_t aTransSiO2 = TMath::Exp(-0.5/AliHMPIDParam::LAbsWin(photEn));
464
465     // evaluate 80 cm of thickness Gap (low density CH4) transparency
466     Double_t aTransGap  = TMath::Exp(-80./AliHMPIDParam::LAbsGap(photEn));
467
468     // evaluate CsI quantum efficiency
469     Double_t aCsIQE            = AliHMPIDParam::QEffCSI(photEn);
470
471     // evaluate total convolution of all material optical properties
472     Double_t aTotConvolution   = aTransRad*aTransSiO2*aTransGap*aCsIQE;
473
474     sEnergProb+=aTotConvolution*photEn;  
475   
476     sProb+=aTotConvolution;  
477 }
478   if(sProb>0) {
479     eMean = sEnergProb/sProb;
480   } else {
481     return DefaultEMean();
482   }
483   Log(Form(" Mean energy photon calculated ---> %f eV ",eMean));
484
485   if(eMean<AliHMPIDParam::EPhotMin() || eMean>AliHMPIDParam::EPhotMax()) return DefaultEMean();
486   
487   return eMean;
488 }   
489 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
490 Double_t AliHMPIDPreprocessor::DefaultEMean()
491 {
492     Double_t eMean = 6.675;      //just set a refractive index for C6F14 at ephot=6.675 eV @ T=25 C
493     AliWarning(Form("Mean energy for photons out of range [%f,%f] in Preprocessor. Default value Eph=%f eV taken.",AliHMPIDParam::EPhotMin(),
494                                                                                                                    AliHMPIDParam::EPhotMax(),
495                                                                                                                    eMean));
496     return eMean;
497 }