]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ANALYSIS/TenderSupplies/AliPHOSTenderSupply.cxx
TOF recalculation with HG only added
[u/mrichter/AliRoot.git] / ANALYSIS / TenderSupplies / AliPHOSTenderSupply.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
17 ///////////////////////////////////////////////////////////////////////////////
18 //                                                                           //
19 //  PHOS tender, recalibrate PHOS clusters                                   //
20 //  and do track matching                                                    //
21 //  Author : Dmitri Peressounko (RRC KI)                                     //
22 //                                                                           //
23 ///////////////////////////////////////////////////////////////////////////////
24
25 #include "TROOT.h"
26 #include "TH2.h"
27 #include "TFile.h"
28
29 #include <AliLog.h>
30 #include <AliVEvent.h>
31 #include <AliAODEvent.h>
32 #include <AliESDEvent.h>
33 #include <AliAnalysisManager.h>
34 #include <AliTender.h>
35 #include <AliCDBManager.h>
36 #include "AliMagF.h"
37 #include "TGeoGlobalMagField.h"
38
39 #include "AliVCluster.h"
40 #include "AliPHOSTenderSupply.h"
41 #include "AliPHOSCalibData.h"
42 #include "AliPHOSGeometry.h"
43 #include "AliPHOSEsdCluster.h"
44 #include "AliPHOSAodCluster.h"
45 #include "AliOADBContainer.h"
46 #include "AliAODCaloCells.h"
47 #include "AliESDCaloCells.h"
48
49 ClassImp(AliPHOSTenderSupply)
50
51 AliPHOSTenderSupply::AliPHOSTenderSupply() :
52   AliTenderSupply()
53   ,fOCDBpass("local://OCDB")
54   ,fNonlinearityVersion("Default")
55   ,fPHOSGeo(0x0)
56   ,fRecoPass(-1)  //to be defined
57   ,fUsePrivateBadMap(0)
58   ,fUsePrivateCalib(0)
59   ,fPHOSCalibData(0x0)
60   ,fTask(0x0)
61 {
62         //
63         // default ctor
64         //
65    for(Int_t i=0;i<10;i++)fNonlinearityParams[i]=0. ;
66    for(Int_t mod=0;mod<5;mod++)fPHOSBadMap[mod]=0x0 ;
67 }
68
69 //_____________________________________________________
70 AliPHOSTenderSupply::AliPHOSTenderSupply(const char *name, const AliTender *tender) :
71   AliTenderSupply(name,tender)
72   ,fOCDBpass("alien:///alice/cern.ch/user/p/prsnko/PHOSrecalibrations/")
73   ,fNonlinearityVersion("Default")
74   ,fPHOSGeo(0x0)
75   ,fRecoPass(-1)  //to be defined
76   ,fUsePrivateBadMap(0)
77   ,fUsePrivateCalib(0)
78   ,fPHOSCalibData(0x0)
79   ,fTask(0x0)
80 {
81         //
82         // named ctor
83         //
84    for(Int_t i=0;i<10;i++)fNonlinearityParams[i]=0. ;
85    for(Int_t mod=0;mod<5;mod++)fPHOSBadMap[mod]=0x0 ;
86 }
87
88 //_____________________________________________________
89 AliPHOSTenderSupply::~AliPHOSTenderSupply()
90 {
91   //Destructor
92   if(fPHOSCalibData)
93     delete fPHOSCalibData;
94   fPHOSCalibData=0x0 ;
95 }
96
97 //_____________________________________________________
98 void AliPHOSTenderSupply::InitTender()
99 {
100   //
101   // Initialise PHOS tender
102   //
103   Int_t runNumber = 0;
104   if(fTender)
105     runNumber = fTender->GetRun();
106   else{
107     if(!fTask){
108       AliError("Neither Tender not Taks was not set") ;
109       return ;
110     }
111     AliAODEvent *aod = dynamic_cast<AliAODEvent*>(fTask->InputEvent()) ;
112     if(aod)
113       runNumber = aod->GetRunNumber() ;
114     else{
115       AliESDEvent *esd = dynamic_cast<AliESDEvent*>(fTask->InputEvent()) ;
116       if(esd)
117         runNumber = esd->GetRunNumber() ;
118       else{
119         AliError("Taks does not contain neither ESD nor AOD") ;
120         return ;
121       }
122     }   
123   }
124     
125   if(fRecoPass<0){ //not defined yet
126     // read if from filename.
127     AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
128     TTree * t = mgr->GetTree();
129     if(t){  
130       TFile * f = t->GetCurrentFile() ;
131       if(f){  
132         TString fname(f->GetName());
133         if(fname.Contains("pass1"))
134            fRecoPass=1;
135         else 
136           if(fname.Contains("pass2"))
137            fRecoPass=2;
138           else 
139             if(fname.Contains("pass3")) 
140               fRecoPass=3;
141       }
142     }
143     if(fRecoPass<0){
144       AliError("Can not find pass number from file name, set it manually");
145     }
146   }
147    
148   //Init geometry 
149   if(!fPHOSGeo){
150     AliOADBContainer geomContainer("phosGeo");
151     geomContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSGeometry.root","PHOSRotationMatrixes");
152     TObjArray *matrixes = (TObjArray*)geomContainer.GetObject(runNumber,"PHOSRotationMatrixes");
153     fPHOSGeo =  AliPHOSGeometry::GetInstance("IHEP") ;
154     for(Int_t mod=0; mod<5; mod++) {
155       if(!matrixes->At(mod)) continue;
156       fPHOSGeo->SetMisalMatrix(((TGeoHMatrix*)matrixes->At(mod)),mod) ;
157       printf(".........Adding Matrix(%d), geo=%p\n",mod,fPHOSGeo) ;
158     }
159   }
160   
161   //Init Bad channels map
162   if(!fUsePrivateBadMap){
163    AliOADBContainer badmapContainer(Form("phosBadMap"));
164     badmapContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSBadMaps.root","phosBadMap");
165     TObjArray *maps = (TObjArray*)badmapContainer.GetObject(runNumber,"phosBadMap");
166     if(!maps){
167       AliError(Form("Can not read Bad map for run %d. \n You may choose to use your map with ForceUsingBadMap()\n",runNumber)) ;    
168     }
169     else{
170       AliInfo(Form("Setting PHOS bad map with name %s \n",maps->GetName())) ;
171       for(Int_t mod=0; mod<5;mod++){
172         if(fPHOSBadMap[mod]) 
173           delete fPHOSBadMap[mod] ;
174         TH2I * h = (TH2I*)maps->At(mod) ;      
175         if(h)
176           fPHOSBadMap[mod]=new TH2I(*h) ;
177       }
178     }    
179   } 
180
181   if(!fUsePrivateCalib){
182     //Init recalibration
183     //Check the pass1-pass2-pass3 reconstruction
184     AliOADBContainer calibContainer("phosRecalibration");
185     calibContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSCalibrations.root","phosRecalibration");
186     TObjArray *recalib = (TObjArray*)calibContainer.GetObject(runNumber,"PHOSRecalibration");
187     if(!recalib){
188       AliFatal(Form("Can not read calibrations for run %d\n. You may choose your specific calibration with ForceUsingCalibration()\n",runNumber)) ;
189     }
190     else{
191       fPHOSCalibData = (AliPHOSCalibData*)recalib->At(fRecoPass-1) ;
192       if(!fPHOSCalibData) {
193         AliFatal(Form("Can not find calibration for run %d, pass %d \n",runNumber, fRecoPass)) ;
194       }
195     }
196   }
197   
198 }
199
200 //_____________________________________________________
201 void AliPHOSTenderSupply::ProcessEvent()
202 {
203   //Choose PHOS clusters and recalibrate them
204   //that it recalculate energy, position and distance 
205   //to closest track extrapolation      
206
207   AliESDEvent *esd = 0x0 ; 
208   AliAODEvent *aod = 0x0 ;
209   if(fTender){
210     esd = fTender->GetEvent();
211     if(!esd)
212       return ;
213   }
214   else{
215     if(!fTask){
216       return ;
217     }
218     esd = dynamic_cast<AliESDEvent*>(fTask->InputEvent()) ;
219     aod = dynamic_cast<AliAODEvent*>(fTask->InputEvent()) ;
220     if(!esd && !aod)
221       return ;
222   }
223     
224   if(!fPHOSCalibData 
225     || (fTender && fTender->RunChanged())){ //In case of Task init called automatically
226     InitTender();
227     
228   }
229
230   TVector3 vertex ;
231   if(esd){
232     const AliESDVertex *esdVertex = esd->GetPrimaryVertex();
233     vertex.SetXYZ(esdVertex->GetX(),esdVertex->GetY(),esdVertex->GetZ());
234   }
235   else{//AOD
236     const AliAODVertex *aodVertex = aod->GetPrimaryVertex();
237     vertex.SetXYZ(aodVertex->GetX(),aodVertex->GetY(),aodVertex->GetZ());
238   }
239   if(vertex.Mag()>99.) //vertex not defined?
240     vertex.SetXYZ(0.,0.,0.) ;
241
242
243   //For re-calibration
244   const Double_t logWeight=4.5 ;  
245
246   if(esd){ //To avoid multiple if in loops we made 
247            //almost identical pecies of code. Please apply changes to both!!!
248     Int_t multClust=esd->GetNumberOfCaloClusters();
249     AliESDCaloCells * cells = esd->GetPHOSCells() ;
250  
251     for (Int_t i=0; i<multClust; i++) {
252       AliESDCaloCluster *clu = esd->GetCaloCluster(i);    
253       if ( !clu->IsPHOS()) continue;
254       
255       Float_t  position[3];
256       clu->GetPosition(position);
257       TVector3 global(position) ;
258       Int_t relId[4] ;
259       fPHOSGeo->GlobalPos2RelId(global,relId) ;
260       Int_t mod  = relId[0] ;
261       Int_t cellX = relId[2];
262       Int_t cellZ = relId[3] ;
263       if ( !IsGoodChannel(mod,cellX,cellZ) ) {
264         clu->SetE(0.) ;
265         continue ;
266       }  
267     
268       //Apply re-Calibreation
269       AliPHOSEsdCluster cluPHOS(*clu);
270       cluPHOS.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
271       cluPHOS.EvalAll(logWeight,vertex);         // recalculate the cluster parameters
272       cluPHOS.SetE(CorrectNonlinearity(cluPHOS.E()));// Users's nonlinearity
273
274       //Correct Misalignment
275       CorrectPHOSMisalignment(global,mod) ;
276       position[0]=global.X() ;
277       position[1]=global.Y() ;
278       position[2]=global.Z() ;
279       cluPHOS.SetPosition(position) ; 
280
281       //Eval CoreDispersion
282       Double_t m02=0.,m20=0.;
283       EvalLambdas(&cluPHOS,m02, m20);   
284       clu->SetDispersion(TestLambda(clu->E(),m20,m02)) ;
285       
286       Float_t  xyz[3];
287       cluPHOS.GetPosition(xyz);
288       clu->SetPosition(xyz);                       //rec.point position in MARS
289       clu->SetE(cluPHOS.E());                           //total or core particle energy
290       //    clu->SetDispersion(cluPHOS.GetDispersion());  //cluster dispersion
291       //    ec->SetPID(rp->GetPID()) ;            //array of particle identification
292       clu->SetM02(cluPHOS.GetM02()) ;               //second moment M2x
293       clu->SetM20(cluPHOS.GetM20()) ;               //second moment M2z
294       Double_t dx=0.,dz=0. ;
295       fPHOSGeo->GlobalPos2RelId(global,relId) ;
296       TVector3 locPos;
297       fPHOSGeo->Global2Local(locPos,global,mod) ;
298
299       Double_t pttrack=0.;
300       Int_t charge=0;
301       FindTrackMatching(mod,&locPos,dx,dz,pttrack,charge) ;
302       Double_t r=TestCPV(dx, dz, pttrack,charge) ;
303       clu->SetTrackDistance(dx,dz); 
304      
305       clu->SetEmcCpvDistance(r);    
306       clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02()));                     //not yet implemented
307       clu->SetTOF(EvalTOFesd(&cluPHOS,cells));       
308
309     }
310   }
311   else{//AOD
312     Int_t multClust=aod->GetNumberOfCaloClusters();
313     AliAODCaloCells * cells = aod->GetPHOSCells() ;
314   
315     for (Int_t i=0; i<multClust; i++) {
316       AliAODCaloCluster *clu = aod->GetCaloCluster(i);    
317       if ( !clu->IsPHOS()) continue;
318       
319       Float_t  position[3];
320       clu->GetPosition(position);
321       TVector3 global(position) ;
322       Int_t relId[4] ;
323       fPHOSGeo->GlobalPos2RelId(global,relId) ;
324       Int_t mod  = relId[0] ;
325       Int_t cellX = relId[2];
326       Int_t cellZ = relId[3] ;
327       if ( !IsGoodChannel(mod,cellX,cellZ) ) {
328         clu->SetE(0.) ;
329         continue ;
330       }  
331       TVector3 locPosOld; //Use it to re-calculate distance to track
332       fPHOSGeo->Global2Local(locPosOld,global,mod) ;
333     
334       //Apply re-Calibreation
335       AliPHOSAodCluster cluPHOS(*clu);
336       cluPHOS.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
337       cluPHOS.EvalAll(logWeight,vertex);         // recalculate the cluster parameters
338       cluPHOS.SetE(CorrectNonlinearity(cluPHOS.E()));// Users's nonlinearity
339
340       //Correct Misalignment
341       cluPHOS.GetPosition(position);
342       global.SetXYZ(position[0],position[1],position[2]);
343       CorrectPHOSMisalignment(global,mod) ;
344       position[0]=global.X() ;
345       position[1]=global.Y() ;
346       position[2]=global.Z() ;
347
348       clu->SetPosition(position);                       //rec.point position in MARS
349       clu->SetE(cluPHOS.E());                           //total or core particle energy
350       clu->SetDispersion(cluPHOS.GetDispersion());  //cluster dispersion
351       //    ec->SetPID(rp->GetPID()) ;            //array of particle identification
352       clu->SetM02(cluPHOS.GetM02()) ;               //second moment M2x
353       clu->SetM20(cluPHOS.GetM20()) ;               //second moment M2z
354       //correct distance to track
355       Double_t dx=clu->GetTrackDx() ;
356       Double_t dz=clu->GetTrackDz() ;
357       if(dx!=-999.){ //there is matched track
358         TVector3 locPos;
359         fPHOSGeo->Global2Local(locPos,global,mod) ;
360         dx+=locPos.X()-locPosOld.X() ;
361         dz+=locPos.Z()-locPosOld.Z() ;      
362         clu->SetTrackDistance(dx,dz);
363       }
364       Double_t r = 999. ; //Big distance
365       int nTracksMatched = clu->GetNTracksMatched();
366       if(nTracksMatched > 0) {
367         AliVTrack* track = dynamic_cast<AliVTrack*> (clu->GetTrackMatched(0));
368         if ( track ) {
369           Double_t pttrack = track->Pt();
370           Short_t charge = track->Charge();
371           r=TestCPV(dx, dz, pttrack,charge) ;
372         }
373       }
374       clu->SetEmcCpvDistance(r); //Distance in sigmas
375      
376       clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02()));                     //not yet implemented
377       clu->SetTOF(EvalTOFaod(&cluPHOS,cells));       
378     }
379   }
380
381 }
382 //___________________________________________________________________________________________________
383 void AliPHOSTenderSupply::FindTrackMatching(Int_t mod,TVector3 *locpos,
384                                             Double_t &dx, Double_t &dz,
385                                             Double_t &pt,Int_t &charge){
386   //Find track with closest extrapolation to cluster
387   AliESDEvent *esd = 0x0 ;
388   if(fTender)
389     esd= fTender->GetEvent();
390   else{ 
391     esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
392   }
393   
394   if(!esd){
395     AliError("ESD is not found") ;
396     return ;
397   }
398   Double_t  magF = esd->GetMagneticField();
399  
400   Double_t magSign = 1.0;
401   if(magF<0)magSign = -1.0;
402   
403   if (!TGeoGlobalMagField::Instance()->GetField()) {
404     AliMagF* field = new AliMagF("Maps","Maps", magSign, magSign, AliMagF::k5kG);
405     TGeoGlobalMagField::Instance()->SetField(field);
406   }
407
408   // *** Start the matching
409   Int_t nt=0;
410   nt = esd->GetNumberOfTracks();
411   //Calculate actual distance to PHOS module
412   TVector3 globaPos ;
413   fPHOSGeo->Local2Global(mod, 0.,0., globaPos) ;
414   const Double_t rPHOS = globaPos.Pt() ; //Distance to center of  PHOS module
415   const Double_t kYmax = 72.+10. ; //Size of the module (with some reserve) in phi direction
416   const Double_t kZmax = 64.+10. ; //Size of the module (with some reserve) in z direction
417   const Double_t kAlpha0=330./180.*TMath::Pi() ; //First PHOS module angular direction
418   const Double_t kAlpha= 20./180.*TMath::Pi() ; //PHOS module angular size
419   Double_t minDistance = 1.e6;
420
421
422   Double_t gposTrack[3] ; 
423
424   Double_t bz = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField();
425   bz = TMath::Sign(0.5*kAlmost0Field,bz) + bz;
426
427   Double_t b[3]; 
428
429     for (Int_t i=0; i<nt; i++) {
430       AliESDtrack *esdTrack=esd->GetTrack(i);
431
432       // Skip the tracks having "wrong" status (has to be checked/tuned)
433       ULong_t status = esdTrack->GetStatus();
434       if ((status & AliESDtrack::kTPCout)   == 0) continue;
435      
436       //Continue extrapolation from TPC outer surface
437       const AliExternalTrackParam *outerParam=esdTrack->GetOuterParam();
438       if (!outerParam) continue;
439       AliExternalTrackParam t(*outerParam);
440      
441       t.GetBxByBz(b) ;
442       //Direction to the current PHOS module
443       Double_t phiMod=kAlpha0-kAlpha*mod ;
444       if(!t.Rotate(phiMod))
445         continue ;
446     
447       Double_t y;                       // Some tracks do not reach the PHOS
448       if (!t.GetYAt(rPHOS,bz,y)) continue; //    because of the bending
449       
450       Double_t z; 
451       if(!t.GetZAt(rPHOS,bz,z))
452         continue ;
453       if (TMath::Abs(z) > kZmax) 
454         continue; // Some tracks miss the PHOS in Z
455       if(TMath::Abs(y) < kYmax){
456         t.PropagateToBxByBz(rPHOS,b);        // Propagate to the matching module
457         //t.CorrectForMaterial(...); // Correct for the TOF material, if needed
458         t.GetXYZ(gposTrack) ;
459         TVector3 globalPositionTr(gposTrack) ;
460         TVector3 localPositionTr ;
461         fPHOSGeo->Global2Local(localPositionTr,globalPositionTr,mod) ;
462         Double_t ddx = locpos->X()-localPositionTr.X();
463         Double_t ddz = locpos->Z()-localPositionTr.Z();
464         Double_t d2 = ddx*ddx + ddz*ddz;
465         if(d2 < minDistance) {
466           dx = ddx ;
467           dz = ddz ;
468           minDistance=d2 ;
469           pt=esdTrack->Pt() ;
470           charge=esdTrack->Charge() ;
471         }
472       }
473     }//Scanned all tracks
474  
475 }
476 //____________________________________________________________
477 Float_t AliPHOSTenderSupply::CorrectNonlinearity(Float_t en){
478
479   //For backward compatibility, if no RecoParameters found
480   if(fNonlinearityVersion=="Default"){
481     return 0.0241+1.0504*en+0.000249*en*en ;
482   }
483
484   if(fNonlinearityVersion=="NoCorrection"){
485     return en ;
486   }
487   if(fNonlinearityVersion=="Gustavo2005"){
488     return fNonlinearityParams[0]+fNonlinearityParams[1]*en + fNonlinearityParams[2]*en*en ;
489   }
490   if(fNonlinearityVersion=="Henrik2010"){
491     return en*(fNonlinearityParams[0]+fNonlinearityParams[1]*TMath::Exp(-en*fNonlinearityParams[2]))*(1.+fNonlinearityParams[3]*TMath::Exp(-en*fNonlinearityParams[4]))*(1.+fNonlinearityParams[6]/(en*en+fNonlinearityParams[5])) ;
492   }
493
494   return en ;
495 }
496 //_____________________________________________________________________________
497 Double_t AliPHOSTenderSupply::TestLambda(Double_t pt,Double_t l1,Double_t l2){
498 //Parameterization for core dispersion   
499 //For R=4.5
500   Double_t   l1Mean  = 1.150200 + 0.097886/(1.+1.486645*pt+0.000038*pt*pt) ;
501   Double_t   l2Mean = 1.574706 + 0.997966*exp(-0.895075*pt)-0.010666*pt ;
502   Double_t   l1Sigma = 0.100255 + 0.337177*exp(-0.517684*pt)+0.001170*pt ;
503   Double_t   l2Sigma = 0.232580 + 0.573401*exp(-0.735903*pt)-0.002325*pt ;
504   Double_t   c = -0.110983 -0.017353/(1.-1.836995*pt+0.934517*pt*pt) ;
505
506 /*
507   //Parameterizatino for full dispersion
508   Double_t l2Mean  = 1.53126+9.50835e+06/(1.+1.08728e+07*pt+1.73420e+06*pt*pt) ;
509   Double_t l1Mean  = 1.12365+0.123770*TMath::Exp(-pt*0.246551)+5.30000e-03*pt ;
510   Double_t l2Sigma = 6.48260e-02+7.60261e+10/(1.+1.53012e+11*pt+5.01265e+05*pt*pt)+9.00000e-03*pt;
511   Double_t l1Sigma = 4.44719e-04+6.99839e-01/(1.+1.22497e+00*pt+6.78604e-07*pt*pt)+9.00000e-03*pt;
512   Double_t c=-0.35-0.550*TMath::Exp(-0.390730*pt) ;
513 */
514   Double_t R2=0.5*(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma + 
515               0.5*(l2-l2Mean)*(l2-l2Mean)/l2Sigma/l2Sigma +
516               0.5*c*(l1-l1Mean)*(l2-l2Mean)/l1Sigma/l2Sigma ;
517   return R2 ;
518   
519 }
520 //____________________________________________________________________________
521 Double_t AliPHOSTenderSupply::TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge){
522   //Parameterization of LHC10h period
523   //_true if neutral_
524   
525   Double_t meanX=0;
526   Double_t meanZ=0.;
527   Double_t sx=TMath::Min(5.4,2.59719e+02*TMath::Exp(-pt/1.02053e-01)+
528               6.58365e-01*5.91917e-01*5.91917e-01/((pt-9.61306e-01)*(pt-9.61306e-01)+5.91917e-01*5.91917e-01)+1.59219);
529   Double_t sz=TMath::Min(2.75,4.90341e+02*1.91456e-02*1.91456e-02/(pt*pt+1.91456e-02*1.91456e-02)+1.60) ;
530   
531   Double_t mf = 0.; //Positive for ++ and negative for --
532   if(fTender){
533     AliESDEvent *esd = fTender->GetEvent();
534     mf = esd->GetMagneticField();
535   }
536   else{ 
537     if(fTask){
538       AliESDEvent *esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
539       if(esd)
540          mf = esd->GetMagneticField();
541       else{
542         AliAODEvent *aod= dynamic_cast<AliAODEvent*>(fTask->InputEvent());
543         if(aod)
544           mf = aod->GetMagneticField();
545       }
546     }else{
547        AliError("Neither Tender nor Task defined") ;    
548     }
549   }
550   
551   if(mf<0.){ //field --
552     meanZ = -0.468318 ;
553     if(charge>0)
554       meanX=TMath::Min(7.3, 3.89994*1.20679*1.20679/(pt*pt+1.20679*1.20679)+0.249029+2.49088e+07*TMath::Exp(-pt*3.33650e+01)) ;
555     else
556       meanX=-TMath::Min(7.7,3.86040*0.912499*0.912499/(pt*pt+0.912499*0.912499)+1.23114+4.48277e+05*TMath::Exp(-pt*2.57070e+01)) ;
557   }
558   else{ //Field ++
559     meanZ= -0.468318;
560     if(charge>0)
561       meanX=-TMath::Min(8.0,3.86040*1.31357*1.31357/(pt*pt+1.31357*1.31357)+0.880579+7.56199e+06*TMath::Exp(-pt*3.08451e+01)) ;
562     else
563       meanX= TMath::Min(6.85, 3.89994*1.16240*1.16240/(pt*pt+1.16240*1.16240)-0.120787+2.20275e+05*TMath::Exp(-pt*2.40913e+01)) ;     
564   }
565
566   Double_t rz=(dz-meanZ)/sz ;
567   Double_t rx=(dx-meanX)/sx ;
568   return TMath::Sqrt(rx*rx+rz*rz) ;
569 }
570
571 //________________________________________________________________________
572 Bool_t AliPHOSTenderSupply::IsGoodChannel(Int_t mod, Int_t ix, Int_t iz)
573 {
574   //Check if this channel belogs to the good ones
575   
576   if(mod>4 || mod<1){
577 //    AliError(Form("No bad map for PHOS module %d ",mod)) ;
578     return kTRUE ;
579   }
580   if(!fPHOSBadMap[mod]){
581 //     AliError(Form("No Bad map for PHOS module %d",mod)) ;
582      return kTRUE ;
583   }
584   if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0)
585     return kFALSE ;
586   else
587     return kTRUE ;
588 }
589 //________________________________________________________________________
590 void AliPHOSTenderSupply::ForceUsingBadMap(const char * filename){
591   //Read TH2I histograms with bad maps from local or alien file 
592   TFile * fbm = TFile::Open(filename) ;
593   if(!fbm || !fbm->IsOpen()){
594     AliError(Form("Can not open BadMaps file %s",filename)) ;
595     return ;
596   }
597   gROOT->cd() ;
598   char key[55] ;
599   for(Int_t mod=1;mod<4; mod++){
600     snprintf(key,55,"PHOS_BadMap_mod%d",mod) ;
601     TH2I * h = (TH2I*)fbm->Get(key) ;
602     if(h)
603        fPHOSBadMap[mod] = new TH2I(*h) ;
604   }
605   fbm->Close() ;
606   fUsePrivateBadMap=kTRUE ;
607 }
608 //________________________________________________________________________
609 void AliPHOSTenderSupply::ForceUsingCalibration(const char * filename){
610   //Read PHOS recalibration parameters from the file.
611   //We assume that file contains single entry: AliPHOSCalibData
612   TFile * fc = TFile::Open(filename) ;
613   if(!fc || !fc->IsOpen()){
614     AliFatal(Form("Can not open Calibration file %s",filename)) ;
615     return ;
616   }
617   fPHOSCalibData = (AliPHOSCalibData*)fc->Get("PHOSCalibration") ;
618   fc->Close() ;
619   fUsePrivateCalib=kTRUE; 
620 }
621 //________________________________________________________________________
622 void AliPHOSTenderSupply::CorrectPHOSMisalignment(TVector3 &global,Int_t mod){
623    //Correct for PHOS modules misalignment 
624   
625     //correct misalignment
626     const Float_t shiftX[6]={0.,-2.3,-2.11,-1.53,0.,0.} ;
627     const Float_t shiftZ[6]={0.,-0.4, 0.52, 0.8,0.,0.} ;
628     TVector3 localPos ;
629     fPHOSGeo->Global2Local(localPos,global,mod) ;
630     fPHOSGeo->Local2Global(mod,localPos.X()+shiftX[mod],localPos.Z()+shiftZ[mod],global);  
631 }
632 //________________________________________________________________________
633 void AliPHOSTenderSupply::EvalLambdas(AliVCluster * clu, Double_t &m02, Double_t &m20){ 
634   //calculate dispecrsion of the cluster in the circle with radius distanceCut around the maximum
635     
636   const Double_t rCut=4.5 ;
637   
638   Double32_t * elist = clu->GetCellsAmplitudeFraction() ;  
639 // Calculates the center of gravity in the local PHOS-module coordinates
640   Float_t wtot = 0;
641   Double_t xc[100]={0} ;
642   Double_t zc[100]={0} ;
643   Double_t x = 0 ;
644   Double_t z = 0 ;
645   Int_t mulDigit=TMath::Min(100,clu->GetNCells()) ;
646   const Double_t logWeight=4.5 ;
647   for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
648     Int_t relid[4] ;
649     Float_t xi ;
650     Float_t zi ;
651     fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
652     fPHOSGeo->RelPosInModule(relid, xi, zi);
653     xc[iDigit]=xi ;
654     zc[iDigit]=zi ;
655     if (clu->E()>0 && elist[iDigit]>0) {
656       Float_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
657       x    += xc[iDigit] * w ;
658       z    += zc[iDigit] * w ;
659       wtot += w ;
660     }
661   }
662   if (wtot>0) {
663     x /= wtot ;
664     z /= wtot ;
665   }
666      
667   wtot = 0. ;
668   Double_t dxx  = 0.;
669   Double_t dzz  = 0.;
670   Double_t dxz  = 0.;
671   Double_t xCut = 0. ;
672   Double_t zCut = 0. ;
673   for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
674     if (clu->E()>0 && elist[iDigit]>0.) {
675         Double_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
676         Double_t xi= xc[iDigit] ;
677         Double_t zi= zc[iDigit] ;
678         if((xi-x)*(xi-x)+(zi-z)*(zi-z) < rCut*rCut){
679           xCut += w * xi ;
680           zCut += w * zi ; 
681           dxx  += w * xi * xi ;
682           dzz  += w * zi * zi ;
683           dxz  += w * xi * zi ; 
684           wtot += w ;
685         }
686     }
687     
688   }
689   if (wtot>0) {
690     xCut/= wtot ;
691     zCut/= wtot ;
692     dxx /= wtot ;
693     dzz /= wtot ;
694     dxz /= wtot ;
695     dxx -= xCut * xCut ;
696     dzz -= zCut * zCut ;
697     dxz -= xCut * zCut ;
698
699     m02 =  0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
700     m20 =  0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz )  ;
701   }
702   else {
703     m20=m02=0.;
704   }
705
706 }
707 //________________________________________________________________________
708 Double_t AliPHOSTenderSupply::EvalTOFesd(AliPHOSEsdCluster * clu,AliESDCaloCells * cells){ 
709   //calculate dispecrsion of the cluster in the circle with radius distanceCut around the maximum
710
711   Double32_t * elist = clu->GetCellsAmplitudeFraction() ;  
712   Int_t mulDigit=clu->GetNCells() ;
713
714   Float_t tMax= 0.; //Time at the maximum
715   Float_t eMax=0. ;
716   for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
717     Int_t absId=clu->GetCellAbsId(iDigit) ;
718     if(cells->GetCellMCLabel(absId)==-2) //This is LG digit. No statistics to calibrate LG timing, remove them from TOF calculation
719       continue ;
720     if( elist[iDigit]>eMax){
721       tMax=CalibrateTOF(cells->GetCellTime(absId),absId,kTRUE) ;
722       eMax=elist[iDigit] ;
723     }
724   }
725  
726   //Try to improve accuracy 
727   //Do not account time of soft cells:
728 //  const Double_t part=0.5 ;
729   Double_t eMin=TMath::Max(0.5,0.5*eMax) ;
730   if(eMin>0.9*eMax)eMin=0.9*eMax ; //At least eMax should remain
731   Float_t wtot = 0.;
732   Double_t t = 0. ;
733   for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
734     Int_t absId=clu->GetCellAbsId(iDigit) ;
735     if(cells->GetCellMCLabel(absId)==-2) //This is LG digit. No statistics to calibrate LG timing, remove them from TOF calculation
736       continue ;
737     
738     Double_t ti=CalibrateTOF(cells->GetCellTime(absId),absId,kTRUE) ;
739     if(TMath::Abs(ti-tMax)>25.e-8) //remove soft cells with wrong time
740       continue ;
741     
742     //Remove too soft cells
743     if(elist[iDigit]<eMin)
744       continue ;
745     
746     if(elist[iDigit]>0){ 
747       //weight = 1./sigma^2
748       //Sigma is parameterization of TOF resolution 16.05.2013
749       Double_t wi2=1./(2.4e-9 + 3.9e-9/elist[iDigit]) ;
750       t+=ti*wi2 ;
751       wtot+=wi2 ;
752     }
753   }
754   if(wtot>0){
755     t=t/wtot ;
756   }  
757   
758   return t ;
759      
760
761 //________________________________________________________________________
762 Double_t AliPHOSTenderSupply::EvalTOFaod(AliPHOSAodCluster * clu,AliAODCaloCells * cells){ 
763   //calculate dispecrsion of the cluster in the circle with radius distanceCut around the maximum
764       
765   Double32_t * elist = clu->GetCellsAmplitudeFraction() ;  
766   Float_t wtot = 0.;
767   Double_t t = 0. ;
768   Int_t mulDigit=clu->GetNCells() ;
769   for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
770     Int_t absId=clu->GetCellAbsId(iDigit) ;
771     Bool_t isHG=kTRUE ;
772     if(cells->GetCellMCLabel(absId)==-2) //This is LG digit. No statistics to calibrate LG timing, remove them from TOF calculation
773       continue ;
774     if(elist[iDigit]>0){ 
775       //weight = 1./sigma^2
776       //Sigma is parameterization of TOF resolution 16.05.2013
777       Double_t wi2=1./(2.4e-9 + 3.9e-9/elist[iDigit]) ;
778       Double_t ti=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
779       t+=ti*wi2 ;
780       wtot+=wi2 ;
781     }
782   }
783   if(wtot>0){
784     t=t/wtot ;
785   }
786   return t ;
787
788 //________________________________________________________________________
789 Double_t AliPHOSTenderSupply::CalibrateTOF(Double_t tof, Int_t absId, Bool_t isHG){
790   //Apply time re-calibration separately for HG and LG channels
791   //By default (if not filled) shifts are zero.  
792     
793   Int_t relId[4];
794   fPHOSGeo->AbsToRelNumbering(absId,relId) ;
795   Int_t   module = relId[0];
796   Int_t   column = relId[3];
797   Int_t   row    = relId[2];
798   if(isHG)
799     tof-=fPHOSCalibData->GetTimeShiftEmc(module, column, row);
800   else
801     tof-=fPHOSCalibData->GetLGTimeShiftEmc(module, column, row);
802  
803   return tof ;
804   
805 }