]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSReconstructor.cxx
Correction (Raphael)
[u/mrichter/AliRoot.git] / PHOS / AliPHOSReconstructor.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 /* $Id$ */
17
18 //_________________________________________________________________________
19 //--
20 //-- Yves Schutz (SUBATECH) 
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and 
22 // derived from STEER/AliReconstructor. 
23 // 
24 // --- ROOT system ---
25 #include "TGeoManager.h"
26 #include "TGeoMatrix.h"
27
28 // --- Standard library ---
29
30 // --- AliRoot header files ---
31 #include "AliLog.h"
32 #include "AliAltroMapping.h"
33 #include "AliESDEvent.h"
34 #include "AliESDCaloCluster.h"
35 #include "AliESDCaloCells.h"
36 #include "AliPHOSReconstructor.h"
37 #include "AliPHOSClusterizerv1.h"
38 #include "AliPHOSTrackSegmentMakerv1.h"
39 #include "AliPHOSPIDv1.h"
40 #include "AliPHOSTracker.h"
41 #include "AliRawReader.h"
42 #include "AliPHOSCalibData.h"
43 #include "AliCDBEntry.h"
44 #include "AliCDBManager.h"
45 #include "AliPHOSTrigger.h"
46 #include "AliPHOSGeometry.h"
47 #include "AliPHOSDigit.h"
48 #include "AliPHOSTrackSegment.h"
49 #include "AliPHOSEmcRecPoint.h"
50 #include "AliPHOSRecParticle.h"
51 #include "AliPHOSRawFitterv0.h"
52 #include "AliPHOSRawFitterv1.h"
53 #include "AliPHOSRawFitterv2.h"
54 #include "AliPHOSRawFitterv3.h"
55 #include "AliPHOSRawFitterv4.h"
56 #include "AliPHOSRawDigiProducer.h"
57 #include "AliPHOSPulseGenerator.h"
58
59 ClassImp(AliPHOSReconstructor)
60
61 Bool_t AliPHOSReconstructor::fgDebug = kFALSE ; 
62 TClonesArray*     AliPHOSReconstructor::fgDigitsArray = 0;   // Array of PHOS digits
63 TObjArray*        AliPHOSReconstructor::fgEMCRecPoints = 0;   // Array of EMC rec.points
64 AliPHOSCalibData * AliPHOSReconstructor::fgCalibData  = 0 ;
65
66
67 //____________________________________________________________________________
68 AliPHOSReconstructor::AliPHOSReconstructor() :
69   fGeom(NULL),fClusterizer(NULL),fTSM(NULL),fPID(NULL),fTmpDigLG(NULL)
70 {
71   // ctor
72   fGeom          = AliPHOSGeometry::GetInstance("IHEP","");
73   fClusterizer   = new AliPHOSClusterizerv1      (fGeom);
74   fTSM           = new AliPHOSTrackSegmentMakerv1(fGeom);
75   fPID           = new AliPHOSPIDv1              (fGeom);
76   fTmpDigLG      = new TClonesArray("AliPHOSDigit",100);
77   fgDigitsArray  = new TClonesArray("AliPHOSDigit",100);
78   fgEMCRecPoints = new TObjArray(100) ;
79   if (!fgCalibData)
80     fgCalibData = new AliPHOSCalibData(-1); //use AliCDBManager's run number
81
82   AliInfo(Form("PHOS bad channel map contains %d bad channel(s).\n",
83                fgCalibData->GetNumOfEmcBadChannels()));
84  
85 }
86
87 //____________________________________________________________________________
88 AliPHOSReconstructor::~AliPHOSReconstructor()
89 {
90   // dtor
91   delete fGeom;
92   delete fClusterizer;
93   delete fTSM;
94   delete fPID;
95   delete fTmpDigLG;
96   delete fgDigitsArray;
97   delete fgEMCRecPoints;
98
99
100 //____________________________________________________________________________
101 void AliPHOSReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
102 {
103   // 'single-event' local reco method called by AliReconstruction; 
104   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
105   // segment maker needs access to the AliESDEvent object to retrieve the tracks reconstructed by 
106   // the global tracking.
107
108   fClusterizer->InitParameters();
109   fClusterizer->SetInput(digitsTree);
110   fClusterizer->SetOutput(clustersTree);
111   if ( Debug() ) 
112     fClusterizer->Digits2Clusters("deb all") ; 
113   else 
114     fClusterizer->Digits2Clusters("") ;
115 }
116
117 //____________________________________________________________________________
118 void AliPHOSReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
119                                    AliESDEvent* esd) const
120 {
121   // This method produces PHOS rec-particles,
122   // then it creates AliESDtracks out of them and
123   // write tracks to the ESD
124
125
126   // do current event; the loop over events is done by AliReconstruction::Run()
127   fTSM->SetESD(esd) ; 
128   fTSM->SetInput(clustersTree);
129   if ( Debug() ) 
130     fTSM->Clusters2TrackSegments("deb all") ;
131   else 
132     fTSM->Clusters2TrackSegments("") ;
133   
134   fPID->SetInput(clustersTree, fTSM->GetTrackSegments()) ; 
135   fPID->SetESD(esd) ; 
136   if ( Debug() ) 
137     fPID->TrackSegments2RecParticles("deb all") ;
138   else 
139     fPID->TrackSegments2RecParticles("") ;
140
141   TClonesArray *recParticles  = fPID->GetRecParticles();
142   Int_t nOfRecParticles = recParticles->GetEntriesFast();
143   
144   AliDebug(2,Form("%d rec. particles, option %s",nOfRecParticles,GetOption()));
145   
146   // Read digits array
147
148   TBranch *branch = digitsTree->GetBranch("PHOS");
149   if (!branch) { 
150     AliError("can't get the branch with the PHOS digits !");
151     return;
152   }
153   branch->SetAddress(&fgDigitsArray);
154   branch->GetEntry(0);
155
156   // Get the clusters array
157
158   TBranch *emcbranch = clustersTree->GetBranch("PHOSEmcRP");
159   if (!emcbranch) { 
160     AliError("can't get the branch with the PHOS EMC clusters !");
161     return;
162   }
163
164   emcbranch->SetAddress(&fgEMCRecPoints);
165   emcbranch->GetEntry(0);
166
167 //   //#########Calculate trigger and set trigger info###########
168
169 //   AliPHOSTrigger tr ;
170 //   //   tr.SetPatchSize(1);//create 4x4 patches
171 //   tr.SetSimulation(kFALSE);
172 //   tr.Trigger(fgDigitsArray);
173   
174 //   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
175 //   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
176 //   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
177 //   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
178
179 //   Int_t iSM2x2      = tr.Get2x2SuperModule();
180 //   Int_t iSMnxn      = tr.GetnxnSuperModule();
181 //   Int_t iCrystalPhi2x2 = tr.Get2x2CrystalPhi();
182 //   Int_t iCrystalPhinxn = tr.GetnxnCrystalPhi();
183 //   Int_t iCrystalEta2x2 = tr.Get2x2CrystalEta();
184 //   Int_t iCrystalEtanxn = tr.GetnxnCrystalEta();
185
186 //   AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",  
187 //                 maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iCrystalPhi2x2, iCrystalEta2x2));
188 //   AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",
189 //                 maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iCrystalPhinxn, iCrystalEtanxn));
190
191 //   // Attention! PHOS modules in order to calculate AbsId need to be 1-5 not 0-4 as returns trigger.
192 //   Int_t iRelId2x2 []= {iSM2x2+1,0,iCrystalPhi2x2,iCrystalEta2x2};
193 //   Int_t iAbsId2x2 =-1;
194 //   Int_t iRelIdnxn []= {iSMnxn+1,0,iCrystalPhinxn,iCrystalEtanxn};
195 //   Int_t iAbsIdnxn =-1;
196 //   TVector3    pos2x2(-1,-1,-1);
197 //   TVector3    posnxn(-1,-1,-1);
198 //   fGeom->RelToAbsNumbering(iRelId2x2, iAbsId2x2);
199 //   fGeom->RelToAbsNumbering(iRelIdnxn, iAbsIdnxn);
200 //   fGeom->RelPosInAlice(iAbsId2x2, pos2x2);
201 //   fGeom->RelPosInAlice(iAbsIdnxn, posnxn);
202
203 //   TArrayF triggerPosition(6);
204 //   triggerPosition[0] = pos2x2(0) ;   
205 //   triggerPosition[1] = pos2x2(1) ;   
206 //   triggerPosition[2] = pos2x2(2) ;  
207 //   triggerPosition[3] = posnxn(0) ;   
208 //   triggerPosition[4] = posnxn(1) ;   
209 //   triggerPosition[5] = posnxn(2) ;  
210
211 //   TArrayF triggerAmplitudes(4);
212 //   triggerAmplitudes[0] = maxAmp2x2 ;   
213 //   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
214 //   triggerAmplitudes[2] = maxAmpnxn ;   
215 //   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
216
217 //   //esd->SetPHOSTriggerCells(triggerPosition);
218 //   esd->AddPHOSTriggerPosition(triggerPosition);
219 //   esd->AddPHOSTriggerAmplitudes(triggerAmplitudes);
220   
221
222   //########################################
223   //############# Fill CaloCells ###########
224   //########################################
225
226   Int_t nDigits = fgDigitsArray->GetEntries();
227   Int_t idignew = 0 ;
228   AliDebug(1,Form("%d digits",nDigits));
229
230   const Int_t knEMC = fGeom->GetNModules()*fGeom->GetNPhi()*fGeom->GetNZ();
231   AliESDCaloCells &phsCells = *(esd->GetPHOSCells());
232   phsCells.CreateContainer(nDigits);
233   phsCells.SetType(AliESDCaloCells::kPHOSCell);
234
235   // Add to CaloCells only EMC digits with non-zero energy 
236   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
237     const AliPHOSDigit * dig = (const AliPHOSDigit*)fgDigitsArray->At(idig);
238     if(dig->GetId() <= knEMC && 
239        Calibrate(dig->GetEnergy(),dig->GetId()) > GetRecoParam()->GetEMCMinE() ){
240       phsCells.SetCell(idignew,dig->GetId(), Calibrate(dig->GetEnergy(),dig->GetId()),
241                                              CalibrateT(dig->GetTime(),dig->GetId()));   
242       idignew++;
243     }
244   }
245   phsCells.SetNumberOfCells(idignew);
246   phsCells.Sort();
247
248   //########################################
249   //############## Fill CaloClusters #######
250   //########################################
251
252   for (Int_t recpart = 0 ; recpart < nOfRecParticles ; recpart++) {
253     AliPHOSRecParticle  *rp    = static_cast<AliPHOSRecParticle*>(recParticles->At(recpart));
254     if (Debug()) 
255       rp->Print();
256     // Get track segment and EMC rec.point associated with this rec.particle
257     AliPHOSTrackSegment *ts    = static_cast<AliPHOSTrackSegment *>(fTSM->GetTrackSegments()
258                                                                     ->At(rp->GetPHOSTSIndex()));
259
260     AliPHOSEmcRecPoint  *emcRP = static_cast<AliPHOSEmcRecPoint *>(fgEMCRecPoints->At(ts->GetEmcIndex()));
261     AliESDCaloCluster   *ec    = new AliESDCaloCluster() ; 
262     
263     Float_t xyz[3];
264     for (Int_t ixyz=0; ixyz<3; ixyz++) 
265       xyz[ixyz] = rp->GetPos()[ixyz];
266     
267     AliDebug(2,Form("Global position xyz=(%f,%f,%f)",xyz[0],xyz[1],xyz[2]));
268    
269     // Create cell lists
270
271     Int_t     cellMult   = emcRP->GetDigitsMultiplicity();
272     Int_t    *digitsList = emcRP->GetDigitsList();
273     Float_t  *rpElist    = emcRP->GetEnergiesList() ;
274     UShort_t *absIdList  = new UShort_t[cellMult];
275     Double_t *fracList   = new Double_t[cellMult];
276
277     for (Int_t iCell=0; iCell<cellMult; iCell++) {
278       AliPHOSDigit *digit = static_cast<AliPHOSDigit *>(fgDigitsArray->At(digitsList[iCell]));
279       absIdList[iCell] = (UShort_t)(digit->GetId());
280       if (digit->GetEnergy() > 0)
281         fracList[iCell] = rpElist[iCell]/(Calibrate(digit->GetEnergy(),digit->GetId()));
282       else
283         fracList[iCell] = 0;
284     }
285
286     //Primaries
287     Int_t  primMult  = 0;
288     Int_t *primList =  emcRP->GetPrimaries(primMult);
289
290     Float_t energy=0.;
291     if (GetRecoParam()->EMCEcore2ESD())
292       energy = emcRP->GetCoreEnergy();
293     else
294       energy = rp->Energy();
295     //Apply nonlinearity correction
296     if(GetRecoParam()->GetEMCEnergyCorrectionOn())
297       energy=CorrectNonlinearity(energy) ;
298
299     // fills the ESDCaloCluster
300     ec->SetType(AliVCluster::kPHOSNeutral);
301     ec->SetPosition(xyz);                       //rec.point position in MARS
302     ec->SetE(energy);                           //total or core particle energy
303     ec->SetDispersion(emcRP->GetDispersion());  //cluster dispersion
304     ec->SetPID(rp->GetPID()) ;            //array of particle identification
305     ec->SetM02(emcRP->GetM2x()) ;               //second moment M2x
306     ec->SetM20(emcRP->GetM2z()) ;               //second moment M2z
307     ec->SetNExMax(emcRP->GetNExMax());          //number of local maxima
308     ec->SetEmcCpvDistance(ts->GetCpvDistance("r")); //Only radius, what about separate x,z????
309     ec->SetTrackDistance(ts->GetCpvDistance("x"),ts->GetCpvDistance("z")); 
310     ec->SetChi2(-1);                     //not yet implemented
311     ec->SetTOF(emcRP->GetTime());               //Time of flight - already calibrated in EMCRecPoint
312
313     //Cells contributing to clusters
314     ec->SetNCells(cellMult);
315     ec->SetCellsAbsId(absIdList);
316     ec->SetCellsAmplitudeFraction(fracList);
317
318     //Distance to the nearest bad crystal
319     ec->SetDistanceToBadChannel(emcRP->GetDistanceToBadCrystal()); 
320   
321     //Array of MC indeces
322     TArrayI arrayPrim(primMult,primList);
323     ec->AddLabels(arrayPrim);
324     
325     //Matched ESD track
326     TArrayI arrayTrackMatched(1);
327     arrayTrackMatched[0]= ts->GetTrackIndex();
328     ec->AddTracksMatched(arrayTrackMatched);
329     
330     Int_t index = esd->AddCaloCluster(ec);
331
332     //Set pointer to this cluster in ESD track
333     Int_t nt=esd->GetNumberOfTracks();
334     for (Int_t itr=0; itr<nt; itr++) {
335       AliESDtrack *esdTrack=esd->GetTrack(itr);
336       if(!esdTrack->IsPHOS())
337         continue ;
338       if(esdTrack->GetPHOScluster()==-recpart){ //we store negative cluster number
339         esdTrack->SetPHOScluster(index) ;
340 //no garatie that only one track matched this cluster
341 //      break ;
342       }
343     }
344  
345     delete ec;   
346     delete [] fracList;
347     delete [] absIdList;
348   }
349   fgDigitsArray ->Clear();
350   fgEMCRecPoints->Clear("C");
351   recParticles  ->Clear();
352
353   //Store PHOS misalignment matrixes
354   FillMisalMatrixes(esd) ;
355
356 }
357
358 //____________________________________________________________________________
359 AliTracker* AliPHOSReconstructor::CreateTracker() const
360 {
361   // creates the PHOS tracker
362   return new AliPHOSTracker();
363 }
364
365 //____________________________________________________________________________
366 void  AliPHOSReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
367 {
368   // Converts raw data to
369   // PHOS digits
370   // Works on a single-event basis
371   rawReader->Reset() ; 
372
373   AliPHOSRawFitterv0 * fitter ;
374
375   const TObjArray* maps = AliPHOSRecoParam::GetMappings();
376   if(!maps) AliFatal("Cannot retrieve ALTRO mappings!!");
377
378   AliAltroMapping *mapping[20];
379   for(Int_t i = 0; i < 20; i++) {
380     mapping[i] = (AliAltroMapping*)maps->At(i);
381   }
382
383   if      (strcmp(GetRecoParam()->EMCFitterVersion(),"v0")==0) 
384     fitter=new AliPHOSRawFitterv0();
385   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v1")==0) 
386     fitter=new AliPHOSRawFitterv1();
387   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v2")==0) 
388     fitter=new AliPHOSRawFitterv2();
389   else if (strcmp(GetRecoParam()->EMCFitterVersion(),"v3")==0) 
390     fitter=new AliPHOSRawFitterv3();
391   else
392     fitter=new AliPHOSRawFitterv4();
393
394   fitter->SubtractPedestals(GetRecoParam()->EMCSubtractPedestals());
395   fitter->SetAmpOffset     (GetRecoParam()->GetGlobalAltroOffset());
396   fitter->SetAmpThreshold  (GetRecoParam()->GetGlobalAltroThreshold());
397
398   TClonesArray *digits = new TClonesArray("AliPHOSDigit",1);
399   digits->SetName("DIGITS");
400   Int_t bufsize = 32000;
401   digitsTree->Branch("PHOS", &digits, bufsize);
402
403   AliPHOSRawDigiProducer rdp(rawReader,mapping);
404
405   rdp.SetEmcMinAmp(GetRecoParam()->GetEMCRawDigitThreshold()); // in ADC
406   rdp.SetCpvMinAmp(GetRecoParam()->GetCPVMinE());
407   rdp.SetSampleQualityCut(GetRecoParam()->GetEMCSampleQualityCut());
408   rdp.MakeDigits(digits,fTmpDigLG,fitter);
409
410   delete fitter ;
411
412   if (AliLog::GetGlobalDebugLevel() == 1) {
413     Int_t modMax=-111;
414     Int_t colMax=-111;
415     Int_t rowMax=-111;
416     Float_t eMax=-333;
417     //!!!for debug!!!
418     
419     Int_t relId[4];
420     for(Int_t iDigit=0; iDigit<digits->GetEntries(); iDigit++) {
421       AliPHOSDigit* digit = (AliPHOSDigit*)digits->At(iDigit);
422       if(digit->GetEnergy()>eMax) {
423         fGeom->AbsToRelNumbering(digit->GetId(),relId);
424         eMax=digit->GetEnergy();
425         modMax=relId[0];
426         rowMax=relId[2];
427         colMax=relId[3];
428       }
429     }
430     
431     AliDebug(1,Form("Digit with max. energy:  modMax %d colMax %d rowMax %d  eMax %f\n\n",
432                     modMax,colMax,rowMax,eMax));
433   }
434
435   digitsTree->Fill();
436   digits->Delete();
437   delete digits;
438 }
439 //==================================================================================
440 Float_t AliPHOSReconstructor::Calibrate(Float_t amp, Int_t absId)const{
441   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
442
443   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
444
445   //Determine rel.position of the cell absolute ID
446   Int_t relId[4];
447   geom->AbsToRelNumbering(absId,relId);
448   Int_t module=relId[0];
449   Int_t row   =relId[2];
450   Int_t column=relId[3];
451   if(relId[1]){ //CPV
452     Float_t calibration = fgCalibData->GetADCchannelCpv(module,column,row);
453     return amp*calibration ;
454   }
455   else{ //EMC
456     Float_t calibration = fgCalibData->GetADCchannelEmc(module,column,row);
457     return amp*calibration ;
458   }
459 }
460 //==================================================================================
461 Float_t AliPHOSReconstructor::CalibrateT(Float_t time, Int_t absId)const{
462   // Calibrate EMC digit, i.e. multiply its Amp by a factor read from CDB
463
464   const AliPHOSGeometry *geom = AliPHOSGeometry::GetInstance() ;
465
466   //Determine rel.position of the cell absolute ID
467   Int_t relId[4];
468   geom->AbsToRelNumbering(absId,relId);
469   Int_t module=relId[0];
470   Int_t row   =relId[2];
471   Int_t column=relId[3];
472   if(relId[1]){ //CPV
473     return 0. ;
474   }
475   else{ //EMC
476     time += fgCalibData->GetTimeShiftEmc(module,column,row);
477     return time ;
478   }
479 }
480 //==================================================================================
481 void AliPHOSReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
482   //Store PHOS matrixes in ESD Header
483
484   //Check, if matrixes was already stored
485   for(Int_t mod=0 ;mod<5; mod++){
486     if(esd->GetPHOSMatrix(mod)!=0)
487       return ;
488   }
489
490   //Create and store matrixes
491   if(!gGeoManager){
492     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
493     return ;
494   }
495   //Note, that owner of copied marixes will be header
496   char path[255] ;
497   TGeoHMatrix * m ;
498   for(Int_t mod=0; mod<5; mod++){
499     snprintf(path,255,"/ALIC_1/PHOS_%d",mod+1) ; //In Geometry modules numbered 1,2,.,5
500     if (gGeoManager->cd(path)){
501       m = gGeoManager->GetCurrentMatrix() ;
502       esd->SetPHOSMatrix(new TGeoHMatrix(*m),mod) ;
503     }
504     else{
505       esd->SetPHOSMatrix(NULL,mod) ;
506     }
507   }
508
509 }
510 //==================================================================================
511 Float_t AliPHOSReconstructor::CorrectNonlinearity(Float_t en){
512
513   //For backward compatibility, if no RecoParameters found
514   if(!GetRecoParam()){
515     return 0.0241+1.0504*en+0.000249*en*en ;
516   }
517
518   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"NoCorrection")==0){
519     return en ;
520   }
521   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Gustavo2005")==0){
522     const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
523     return par[0]+par[1]*en + par[2]*en*en ;
524   }
525   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"Henrik2010")==0){
526      const Float_t *par=GetRecoParam()->GetNonlinearityParams() ;
527      return en*(par[0]+par[1]*TMath::Exp(-en*par[2]))*(1.+par[3]*TMath::Exp(-en*par[4]))*(1.+par[6]/(en*en+par[5])) ;
528   }
529   //For backward compatibility
530   if(strcmp(GetRecoParam()->GetNonlinearityCorrectionVersion(),"")==0){
531     return 0.0241+1.0504*en+0.000249*en*en ;
532   }
533   return en ;
534 }
535