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