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