1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
17 ///////////////////////////////////////////////////////////////////////////////
19 // PHOS tender, recalibrate PHOS clusters //
20 // and do track matching //
21 // Author : Dmitri Peressounko (RRC KI) //
23 ///////////////////////////////////////////////////////////////////////////////
30 #include <AliVEvent.h>
31 #include <AliAODEvent.h>
32 #include <AliESDEvent.h>
33 #include <AliAnalysisManager.h>
34 #include <AliTender.h>
35 #include <AliCDBManager.h>
37 #include "TGeoGlobalMagField.h"
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"
49 ClassImp(AliPHOSTenderSupply)
51 AliPHOSTenderSupply::AliPHOSTenderSupply() :
53 ,fOCDBpass("local://OCDB")
54 ,fNonlinearityVersion("Default")
56 ,fRecoPass(-1) //to be defined
67 for(Int_t i=0;i<10;i++)fNonlinearityParams[i]=0. ;
68 for(Int_t mod=0;mod<5;mod++)fPHOSBadMap[mod]=0x0 ;
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")
77 ,fRecoPass(-1) //to be defined
88 for(Int_t i=0;i<10;i++)fNonlinearityParams[i]=0. ;
89 for(Int_t mod=0;mod<5;mod++)fPHOSBadMap[mod]=0x0 ;
92 //_____________________________________________________
93 AliPHOSTenderSupply::~AliPHOSTenderSupply()
97 delete fPHOSCalibData;
101 //_____________________________________________________
102 void AliPHOSTenderSupply::InitTender()
105 // Initialise PHOS tender
109 runNumber = fTender->GetRun();
112 AliError("Neither Tender not Taks was not set") ;
115 AliAODEvent *aod = dynamic_cast<AliAODEvent*>(fTask->InputEvent()) ;
117 runNumber = aod->GetRunNumber() ;
119 AliESDEvent *esd = dynamic_cast<AliESDEvent*>(fTask->InputEvent()) ;
121 runNumber = esd->GetRunNumber() ;
123 AliError("Taks does not contain neither ESD nor AOD") ;
129 //In MC always reco pass 1
133 if(fRecoPass<0){ //not defined yet
134 // read if from filename.
135 AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager();
136 TTree * t = mgr->GetTree();
138 TFile * f = t->GetCurrentFile() ;
140 TString fname(f->GetName());
141 if(fname.Contains("pass1"))
144 if(fname.Contains("pass2"))
147 if(fname.Contains("pass3"))
150 if(fname.Contains("pass4"))
155 AliError("Can not find pass number from file name, set it manually");
161 AliOADBContainer geomContainer("phosGeo");
162 if(fIsMC) //use excatly the same geometry as in simulation
163 geomContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSMCGeometry.root","PHOSMCRotationMatrixes");
164 else //Use best approaximation to real geometry
165 geomContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSGeometry.root","PHOSRotationMatrixes");
166 TObjArray *matrixes = (TObjArray*)geomContainer.GetObject(runNumber,"PHOSRotationMatrixes");
167 fPHOSGeo = AliPHOSGeometry::GetInstance("IHEP") ;
168 for(Int_t mod=0; mod<5; mod++) {
169 if(!matrixes->At(mod)) continue;
170 fPHOSGeo->SetMisalMatrix(((TGeoHMatrix*)matrixes->At(mod)),mod) ;
171 printf(".........Adding Matrix(%d), geo=%p\n",mod,fPHOSGeo) ;
172 ((TGeoHMatrix*)matrixes->At(mod))->Print() ;
176 //Init Bad channels map
177 if(!fUsePrivateBadMap){
178 AliOADBContainer badmapContainer(Form("phosBadMap"));
179 badmapContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSBadMaps.root","phosBadMap");
180 TObjArray *maps = (TObjArray*)badmapContainer.GetObject(runNumber,"phosBadMap");
182 AliError(Form("Can not read Bad map for run %d. \n You may choose to use your map with ForceUsingBadMap()\n",runNumber)) ;
185 AliInfo(Form("Setting PHOS bad map with name %s \n",maps->GetName())) ;
186 for(Int_t mod=0; mod<5;mod++){
188 delete fPHOSBadMap[mod] ;
189 TH2I * h = (TH2I*)maps->At(mod) ;
191 fPHOSBadMap[mod]=new TH2I(*h) ;
196 if(!fUsePrivateCalib){
197 if(fIsMC){ //re/de-calibration for MC productions
199 AliOADBContainer calibContainer("phosRecalibration");
200 calibContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSMCCalibrations.root","phosRecalibration");
202 TObjArray *recalib = (TObjArray*)calibContainer.GetObject(runNumber,"PHOSRecalibration");
204 AliFatal(Form("Can not read calibrations for run %d\n. You may choose your specific calibration with ForceUsingCalibration()\n",runNumber)) ;
207 //Now try to find object with proper name
208 for(Int_t i=0; i<recalib->GetEntriesFast(); i++){
209 AliPHOSCalibData * tmp = (AliPHOSCalibData*)recalib->At(i) ;
210 if(fMCProduction.CompareTo(tmp->GetName())==0){
211 fPHOSCalibData = tmp ;
215 if(!fPHOSCalibData) {
216 AliFatal(Form("Can not find calibration for run %d, and name %s \n",runNumber, fMCProduction.Data())) ;
223 //Check the pass1-pass2-pass3 reconstruction
224 AliOADBContainer calibContainer("phosRecalibration");
225 calibContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSCalibrations.root","phosRecalibration");
226 TObjArray *recalib = (TObjArray*)calibContainer.GetObject(runNumber,"PHOSRecalibration");
228 AliFatal(Form("Can not read calibrations for run %d\n. You may choose your specific calibration with ForceUsingCalibration()\n",runNumber)) ;
231 fPHOSCalibData = (AliPHOSCalibData*)recalib->At(fRecoPass-1) ;
232 if(!fPHOSCalibData) {
233 AliFatal(Form("Can not find calibration for run %d, pass %d \n",runNumber, fRecoPass)) ;
241 //_____________________________________________________
242 void AliPHOSTenderSupply::ProcessEvent()
244 //Choose PHOS clusters and recalibrate them
245 //that it recalculate energy, position and distance
246 //to closest track extrapolation
248 AliESDEvent *esd = 0x0 ;
249 AliAODEvent *aod = 0x0 ;
251 esd = fTender->GetEvent();
259 esd = dynamic_cast<AliESDEvent*>(fTask->InputEvent()) ;
260 aod = dynamic_cast<AliAODEvent*>(fTask->InputEvent()) ;
266 || (fTender && fTender->RunChanged())){ //In case of Task init called automatically
273 const AliESDVertex *esdVertex = esd->GetPrimaryVertex();
274 vertex.SetXYZ(esdVertex->GetX(),esdVertex->GetY(),esdVertex->GetZ());
277 const AliAODVertex *aodVertex = aod->GetPrimaryVertex();
278 vertex.SetXYZ(aodVertex->GetX(),aodVertex->GetY(),aodVertex->GetZ());
280 if(vertex.Mag()>99.) //vertex not defined?
281 vertex.SetXYZ(0.,0.,0.) ;
285 const Double_t logWeight=4.5 ;
287 if(esd){ //To avoid multiple if in loops we made
288 //almost identical pecies of code. Please apply changes to both!!!
289 Int_t multClust=esd->GetNumberOfCaloClusters();
290 AliESDCaloCells * cells = esd->GetPHOSCells() ;
292 for (Int_t i=0; i<multClust; i++) {
293 AliESDCaloCluster *clu = esd->GetCaloCluster(i);
294 if ( !clu->IsPHOS()) continue;
297 //Apply re-Calibreation
298 AliPHOSEsdCluster cluPHOS(*clu);
299 cluPHOS.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
300 cluPHOS.EvalAll(logWeight,vertex); // recalculate the cluster parameters
301 cluPHOS.SetE(CorrectNonlinearity(cluPHOS.E()));// Users's nonlinearity
304 cluPHOS.GetPosition(position);
305 clu->SetPosition(position); //rec.point position in MARS
306 TVector3 global(position) ;
308 fPHOSGeo->GlobalPos2RelId(global,relId) ;
309 Int_t mod = relId[0] ;
310 Int_t cellX = relId[2];
311 Int_t cellZ = relId[3] ;
312 if ( !IsGoodChannel(mod,cellX,cellZ) ) {
317 Double_t ecore=CoreEnergy(&cluPHOS) ;
318 ecore=CorrectNonlinearity(ecore) ;
320 clu->SetE(cluPHOS.E()); //total particle energy
321 clu->SetCoreEnergy(ecore); //core particle energy
323 //Eval FullDispersion
324 clu->SetDispersion(TestFullLambda(clu->E(),cluPHOS.GetM20(),cluPHOS.GetM02())) ;
325 //Eval CoreDispersion
326 Double_t m02=0.,m20=0.;
327 EvalLambdas(&cluPHOS,m02, m20);
328 clu->SetChi2(TestCoreLambda(clu->E(),m20,m02)); //not yet implemented
329 clu->SetM02(m02) ; //second moment M2x
330 clu->SetM20(m20) ; //second moment M2z
332 //correct distance to track
333 Double_t dx=0.,dz=0. ;
334 fPHOSGeo->GlobalPos2RelId(global,relId) ;
336 fPHOSGeo->Global2Local(locPos,global,mod) ;
340 FindTrackMatching(mod,&locPos,dx,dz,pttrack,charge) ;
341 Double_t r=TestCPV(dx, dz, pttrack,charge) ;
342 clu->SetTrackDistance(dx,dz);
344 clu->SetEmcCpvDistance(r);
346 Double_t tof=EvalTOF(&cluPHOS,cells);
347 // if(TMath::Abs(tof-clu->GetTOF())>100.e-9) //something wrong in cell TOF!
348 // tof=clu->GetTOF() ;
350 Double_t minDist=clu->GetDistanceToBadChannel() ;//Already calculated
351 DistanceToBadChannel(mod,&locPos,minDist);
352 clu->SetDistanceToBadChannel(minDist) ;
354 Double_t ecross = EvalEcross(&cluPHOS);
355 clu->SetMCEnergyFraction(ecross) ;
359 Int_t multClust=aod->GetNumberOfCaloClusters();
360 AliAODCaloCells * cells = aod->GetPHOSCells() ;
362 for (Int_t i=0; i<multClust; i++) {
363 AliAODCaloCluster *clu = aod->GetCaloCluster(i);
364 if ( !clu->IsPHOS()) continue;
367 //Apply re-Calibreation
368 AliPHOSAodCluster cluPHOS(*clu);
369 cluPHOS.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
370 cluPHOS.EvalAll(logWeight,vertex); // recalculate the cluster parameters
371 cluPHOS.SetE(CorrectNonlinearity(cluPHOS.E()));// Users's nonlinearity
374 cluPHOS.GetPosition(position);
375 clu->SetPosition(position); //rec.point position in MARS
376 TVector3 global(position) ;
378 fPHOSGeo->GlobalPos2RelId(global,relId) ;
379 Int_t mod = relId[0] ;
380 Int_t cellX = relId[2];
381 Int_t cellZ = relId[3] ;
382 if ( !IsGoodChannel(mod,cellX,cellZ) ) {
386 TVector3 locPosOld; //Use it to re-calculate distance to track
387 fPHOSGeo->Global2Local(locPosOld,global,mod) ;
389 Double_t ecore=CoreEnergy(&cluPHOS) ;
390 ecore=CorrectNonlinearity(ecore) ;
392 clu->SetE(cluPHOS.E()); //total particle energy
393 clu->SetCoreEnergy(ecore); //core particle energy
395 //Eval FullDispersion
396 clu->SetDispersion(TestFullLambda(clu->E(),cluPHOS.GetM20(),cluPHOS.GetM02())) ;
397 //Eval CoreDispersion
398 Double_t m02=0.,m20=0.;
399 EvalLambdas(&cluPHOS,m02, m20);
400 clu->SetChi2(TestCoreLambda(clu->E(),m20,m02)); //not yet implemented
401 clu->SetM02(m02) ; //second moment M2x
402 clu->SetM20(m20) ; //second moment M2z
404 //correct distance to track
405 Double_t dx=clu->GetTrackDx() ;
406 Double_t dz=clu->GetTrackDz() ;
408 fPHOSGeo->Global2Local(locPos,global,mod) ;
409 if(dx!=-999.){ //there is matched track
410 dx+=locPos.X()-locPosOld.X() ;
411 dz+=locPos.Z()-locPosOld.Z() ;
412 clu->SetTrackDistance(dx,dz);
414 Double_t r = 999. ; //Big distance
415 int nTracksMatched = clu->GetNTracksMatched();
416 if(nTracksMatched > 0) {
417 AliVTrack* track = dynamic_cast<AliVTrack*> (clu->GetTrackMatched(0));
419 Double_t pttrack = track->Pt();
420 Short_t charge = track->Charge();
421 r=TestCPV(dx, dz, pttrack,charge) ;
424 clu->SetEmcCpvDistance(r); //Distance in sigmas
427 Double_t tof=EvalTOF(&cluPHOS,cells);
428 // if(TMath::Abs(tof-clu->GetTOF())>100.e-9) //something wrong in cell TOF!
429 // tof=clu->GetTOF() ;
431 Double_t minDist=clu->GetDistanceToBadChannel() ;//Already calculated
432 DistanceToBadChannel(mod,&locPos,minDist);
433 clu->SetDistanceToBadChannel(minDist) ;
435 Double_t ecross = EvalEcross(&cluPHOS);
436 clu->SetMCEnergyFraction(ecross) ;
441 //___________________________________________________________________________________________________
442 void AliPHOSTenderSupply::FindTrackMatching(Int_t mod,TVector3 *locpos,
443 Double_t &dx, Double_t &dz,
444 Double_t &pt,Int_t &charge){
445 //Find track with closest extrapolation to cluster
446 AliESDEvent *esd = 0x0 ;
448 esd= fTender->GetEvent();
450 esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
454 AliError("ESD is not found") ;
457 Double_t magF = esd->GetMagneticField();
459 Double_t magSign = 1.0;
460 if(magF<0)magSign = -1.0;
462 if (!TGeoGlobalMagField::Instance()->GetField()) {
463 AliError("Margnetic filed was not initialized, use default") ;
464 AliMagF* field = new AliMagF("Maps","Maps", magSign, magSign, AliMagF::k5kG);
465 TGeoGlobalMagField::Instance()->SetField(field);
468 // *** Start the matching
470 nt = esd->GetNumberOfTracks();
471 //Calculate actual distance to PHOS module
473 fPHOSGeo->Local2Global(mod, 0.,0., globaPos) ;
474 const Double_t rPHOS = globaPos.Pt() ; //Distance to center of PHOS module
475 const Double_t kYmax = 72.+10. ; //Size of the module (with some reserve) in phi direction
476 const Double_t kZmax = 64.+10. ; //Size of the module (with some reserve) in z direction
477 const Double_t kAlpha0=330./180.*TMath::Pi() ; //First PHOS module angular direction
478 const Double_t kAlpha= 20./180.*TMath::Pi() ; //PHOS module angular size
479 Double_t minDistance = 1.e6;
482 Double_t gposTrack[3] ;
484 Double_t bz = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField();
485 bz = TMath::Sign(0.5*kAlmost0Field,bz) + bz;
489 for (Int_t i=0; i<nt; i++) {
490 AliESDtrack *esdTrack=esd->GetTrack(i);
492 // Skip the tracks having "wrong" status (has to be checked/tuned)
493 ULong_t status = esdTrack->GetStatus();
494 if ((status & AliESDtrack::kTPCout) == 0) continue;
496 //Continue extrapolation from TPC outer surface
497 const AliExternalTrackParam *outerParam=esdTrack->GetOuterParam();
498 if (!outerParam) continue;
499 AliExternalTrackParam t(*outerParam);
502 //Direction to the current PHOS module
503 Double_t phiMod=kAlpha0-kAlpha*mod ;
504 if(!t.Rotate(phiMod))
507 Double_t y; // Some tracks do not reach the PHOS
508 if (!t.GetYAt(rPHOS,bz,y)) continue; // because of the bending
511 if(!t.GetZAt(rPHOS,bz,z))
513 if (TMath::Abs(z) > kZmax)
514 continue; // Some tracks miss the PHOS in Z
515 if(TMath::Abs(y) < kYmax){
516 t.PropagateToBxByBz(rPHOS,b); // Propagate to the matching module
517 //t.CorrectForMaterial(...); // Correct for the TOF material, if needed
518 t.GetXYZ(gposTrack) ;
519 TVector3 globalPositionTr(gposTrack) ;
520 TVector3 localPositionTr ;
521 fPHOSGeo->Global2Local(localPositionTr,globalPositionTr,mod) ;
522 Double_t ddx = locpos->X()-localPositionTr.X();
523 Double_t ddz = locpos->Z()-localPositionTr.Z();
524 Double_t d2 = ddx*ddx + ddz*ddz;
525 if(d2 < minDistance) {
530 charge=esdTrack->Charge() ;
533 }//Scanned all tracks
536 //____________________________________________________________
537 Float_t AliPHOSTenderSupply::CorrectNonlinearity(Float_t en){
539 //For backward compatibility, if no RecoParameters found
540 if(fNonlinearityVersion=="Default"){
541 return 0.0241+1.0504*en+0.000249*en*en ;
543 if(fNonlinearityVersion=="MC"){ //Default + some correction
544 return (0.0241+1.0504*en+0.000249*en*en)*fNonlinearityParams[0]*(1+fNonlinearityParams[1]/(1.+en*en/fNonlinearityParams[2]/fNonlinearityParams[2])) ;
547 if(fNonlinearityVersion=="NoCorrection"){
550 if(fNonlinearityVersion=="Gustavo2005"){
551 return fNonlinearityParams[0]+fNonlinearityParams[1]*en + fNonlinearityParams[2]*en*en ;
553 if(fNonlinearityVersion=="Henrik2010"){
554 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])) ;
559 //_____________________________________________________________________________
560 Double_t AliPHOSTenderSupply::TestCoreLambda(Double_t pt,Double_t l1,Double_t l2){
561 //Parameterization for core dispersion
563 Double_t l1Mean = 1.150200 + 0.097886/(1.+1.486645*pt+0.000038*pt*pt) ;
564 Double_t l2Mean = 1.574706 + 0.997966*exp(-0.895075*pt)-0.010666*pt ;
565 Double_t l1Sigma = 0.100255 + 0.337177*exp(-0.517684*pt)+0.001170*pt ;
566 Double_t l2Sigma = 0.232580 + 0.573401*exp(-0.735903*pt)-0.002325*pt ;
567 Double_t c = -0.110983 -0.017353/(1.-1.836995*pt+0.934517*pt*pt) ;
570 //Parameterizatino for full dispersion
571 Double_t l2Mean = 1.53126+9.50835e+06/(1.+1.08728e+07*pt+1.73420e+06*pt*pt) ;
572 Double_t l1Mean = 1.12365+0.123770*TMath::Exp(-pt*0.246551)+5.30000e-03*pt ;
573 Double_t l2Sigma = 6.48260e-02+7.60261e+10/(1.+1.53012e+11*pt+5.01265e+05*pt*pt)+9.00000e-03*pt;
574 Double_t l1Sigma = 4.44719e-04+6.99839e-01/(1.+1.22497e+00*pt+6.78604e-07*pt*pt)+9.00000e-03*pt;
575 Double_t c=-0.35-0.550*TMath::Exp(-0.390730*pt) ;
577 Double_t R2=0.5*(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma +
578 0.5*(l2-l2Mean)*(l2-l2Mean)/l2Sigma/l2Sigma +
579 0.5*c*(l1-l1Mean)*(l2-l2Mean)/l1Sigma/l2Sigma ;
583 //_____________________________________________________________________________
584 Double_t AliPHOSTenderSupply::TestFullLambda(Double_t pt,Double_t l1,Double_t l2){
585 //Parameterization for full dispersion
586 //Parameterizatino for full dispersion
587 Double_t l2Mean = 1.53126+9.50835e+06/(1.+1.08728e+07*pt+1.73420e+06*pt*pt) ;
588 Double_t l1Mean = 1.12365+0.123770*TMath::Exp(-pt*0.246551)+5.30000e-03*pt ;
589 Double_t l2Sigma = 6.48260e-02+7.60261e+10/(1.+1.53012e+11*pt+5.01265e+05*pt*pt)+9.00000e-03*pt;
590 Double_t l1Sigma = 4.44719e-04+6.99839e-01/(1.+1.22497e+00*pt+6.78604e-07*pt*pt)+9.00000e-03*pt;
591 Double_t c=-0.35-0.550*TMath::Exp(-0.390730*pt) ;
593 Double_t R2=0.5*(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma +
594 0.5*(l2-l2Mean)*(l2-l2Mean)/l2Sigma/l2Sigma +
595 0.5*c*(l1-l1Mean)*(l2-l2Mean)/l1Sigma/l2Sigma ;
599 //____________________________________________________________________________
600 Double_t AliPHOSTenderSupply::TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge){
601 //Parameterization of LHC10h period
606 Double_t sx=TMath::Min(5.4,2.59719e+02*TMath::Exp(-pt/1.02053e-01)+
607 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);
608 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) ;
610 Double_t mf = 0.; //Positive for ++ and negative for --
612 AliESDEvent *esd = fTender->GetEvent();
613 mf = esd->GetMagneticField();
617 AliESDEvent *esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
619 mf = esd->GetMagneticField();
621 AliAODEvent *aod= dynamic_cast<AliAODEvent*>(fTask->InputEvent());
623 mf = aod->GetMagneticField();
626 AliError("Neither Tender nor Task defined") ;
630 if(mf<0.){ //field --
633 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)) ;
635 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)) ;
640 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)) ;
642 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)) ;
645 Double_t rz=(dz-meanZ)/sz ;
646 Double_t rx=(dx-meanX)/sx ;
647 return TMath::Sqrt(rx*rx+rz*rz) ;
650 //________________________________________________________________________
651 Bool_t AliPHOSTenderSupply::IsGoodChannel(Int_t mod, Int_t ix, Int_t iz)
653 //Check if this channel belogs to the good ones
656 // AliError(Form("No bad map for PHOS module %d ",mod)) ;
659 if(!fPHOSBadMap[mod]){
660 // AliError(Form("No Bad map for PHOS module %d",mod)) ;
663 if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0)
668 //________________________________________________________________________
669 void AliPHOSTenderSupply::ForceUsingBadMap(const char * filename){
670 //Read TH2I histograms with bad maps from local or alien file
671 TFile * fbm = TFile::Open(filename) ;
672 if(!fbm || !fbm->IsOpen()){
673 AliError(Form("Can not open BadMaps file %s",filename)) ;
678 for(Int_t mod=1;mod<4; mod++){
679 snprintf(key,55,"PHOS_BadMap_mod%d",mod) ;
680 TH2I * h = (TH2I*)fbm->Get(key) ;
682 fPHOSBadMap[mod] = new TH2I(*h) ;
685 fUsePrivateBadMap=kTRUE ;
687 //________________________________________________________________________
688 void AliPHOSTenderSupply::ForceUsingCalibration(const char * filename){
689 //Read PHOS recalibration parameters from the file.
690 //We assume that file contains single entry: AliPHOSCalibData
691 TFile * fc = TFile::Open(filename) ;
692 if(!fc || !fc->IsOpen()){
693 AliFatal(Form("Can not open Calibration file %s",filename)) ;
696 fPHOSCalibData = (AliPHOSCalibData*)fc->Get("PHOSCalibration") ;
698 fUsePrivateCalib=kTRUE;
700 //________________________________________________________________________
701 void AliPHOSTenderSupply::CorrectPHOSMisalignment(TVector3 &global,Int_t mod){
702 //Correct for PHOS modules misalignment
704 //correct misalignment
705 const Float_t shiftX[6]={0.,-2.3,-2.11,-1.53,0.,0.} ;
706 const Float_t shiftZ[6]={0.,-0.4, 0.52, 0.8,0.,0.} ;
708 fPHOSGeo->Global2Local(localPos,global,mod) ;
709 fPHOSGeo->Local2Global(mod,localPos.X()+shiftX[mod],localPos.Z()+shiftZ[mod],global);
711 //________________________________________________________________________
712 void AliPHOSTenderSupply::EvalLambdas(AliVCluster * clu, Double_t &m02, Double_t &m20){
713 //calculate dispecrsion of the cluster in the circle with radius distanceCut around the maximum
715 const Double_t rCut=4.5 ;
717 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
718 // Calculates the center of gravity in the local PHOS-module coordinates
720 Double_t xc[100]={0} ;
721 Double_t zc[100]={0} ;
724 Int_t mulDigit=TMath::Min(100,clu->GetNCells()) ;
725 const Double_t logWeight=4.5 ;
726 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
730 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
731 fPHOSGeo->RelPosInModule(relid, xi, zi);
734 if (clu->E()>0 && elist[iDigit]>0) {
735 Float_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
736 x += xc[iDigit] * w ;
737 z += zc[iDigit] * w ;
752 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
753 if (clu->E()>0 && elist[iDigit]>0.) {
754 Double_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
755 Double_t xi= xc[iDigit] ;
756 Double_t zi= zc[iDigit] ;
757 if((xi-x)*(xi-x)+(zi-z)*(zi-z) < rCut*rCut){
778 m02 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
779 m20 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
786 //____________________________________________________________________________
787 Double_t AliPHOSTenderSupply::CoreEnergy(AliVCluster * clu){
788 //calculate energy of the cluster in the circle with radius distanceCut around the maximum
790 //Can not use already calculated coordinates?
791 //They have incidence correction...
792 const Double_t distanceCut =3.5 ;
793 const Double_t logWeight=4.5 ;
795 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
796 // Calculates the center of gravity in the local PHOS-module coordinates
798 Double_t xc[100]={0} ;
799 Double_t zc[100]={0} ;
802 Int_t mulDigit=TMath::Min(100,clu->GetNCells()) ;
803 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
807 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
808 fPHOSGeo->RelPosInModule(relid, xi, zi);
811 if (clu->E()>0 && elist[iDigit]>0) {
812 Float_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
813 x += xc[iDigit] * w ;
814 z += zc[iDigit] * w ;
823 for(Int_t iDigit=0; iDigit < mulDigit; iDigit++) {
824 Double_t distance = TMath::Sqrt((xc[iDigit]-x)*(xc[iDigit]-x)+(zc[iDigit]-z)*(zc[iDigit]-z)) ;
825 if(distance < distanceCut)
826 coreE += elist[iDigit] ;
828 //Apply non-linearity correction
831 //____________________________________________________________________________
832 Double_t AliPHOSTenderSupply::EvalEcross(AliVCluster * clu){
833 //Calculate propoerion of the cluster energy in cross around the
834 //cell with maximal energy deposition. Can be used to reject exotic clusters
836 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
837 Int_t mulDigit=clu->GetNCells() ;
838 // Calculates the center of gravity in the local PHOS-module coordinates
839 //Find cell with max E
842 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
843 if(elist[iDigit]>eMax){
848 //Calculate e in cross
851 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iMax), relidMax) ;
853 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
854 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
855 if(TMath::Abs(relid[2]-relidMax[2])+TMath::Abs(relid[3]-relidMax[3])==1)
856 eCross+= elist[iDigit] ;
859 return 1.-eCross/eMax ;
865 //________________________________________________________________________
866 Double_t AliPHOSTenderSupply::EvalTOF(AliVCluster * clu,AliVCaloCells * cells){
867 //Evaluate TOF of the cluster after re-calibration
868 //TOF here is weighted average of digits
869 // -within 50ns from the most energetic cell
872 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
873 Int_t mulDigit=clu->GetNCells() ;
875 Float_t tMax= 0.; //Time at the maximum
877 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
878 Int_t absId=clu->GetCellAbsId(iDigit) ;
879 Bool_t isHG=cells->GetCellHighGain(absId) ;
880 if( elist[iDigit]>eMax){
881 tMax=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
886 //Try to improve accuracy
887 //Do not account time of soft cells:
888 // const Double_t part=0.5 ;
889 Double_t eMin=TMath::Min(0.5,0.2*eMax) ;
892 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
893 Int_t absId=clu->GetCellAbsId(iDigit) ;
894 Bool_t isHG=cells->GetCellHighGain(absId) ;
896 Double_t ti=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
897 if(TMath::Abs(ti-tMax)>50.e-9) //remove soft cells with wrong time
900 //Remove too soft cells
901 if(elist[iDigit]<eMin)
905 //weight = 1./sigma^2
906 //Sigma is parameterization of TOF resolution 16.05.2013
909 wi2=1./(2.4 + 3.9/elist[iDigit]) ;
911 wi2=1./(2.4 + 3.9/(0.1*elist[iDigit])) ; //E of LG digit is 1/16 of correcponding HG
926 //________________________________________________________________________
927 Double_t AliPHOSTenderSupply::CalibrateTOF(Double_t tof, Int_t absId, Bool_t isHG){
928 //Apply time re-calibration separately for HG and LG channels
929 //By default (if not filled) shifts are zero.
932 fPHOSGeo->AbsToRelNumbering(absId,relId) ;
933 Int_t module = relId[0];
934 Int_t column = relId[3];
935 Int_t row = relId[2];
937 tof-=fPHOSCalibData->GetTimeShiftEmc(module, column, row);
939 tof-=fPHOSCalibData->GetLGTimeShiftEmc(module, column, row);
944 //________________________________________________________________________
945 void AliPHOSTenderSupply::DistanceToBadChannel(Int_t mod, TVector3 * locPos, Double_t &minDist){
946 //Check if distance to bad channel was reduced
947 Int_t range = minDist/2.2 +1 ; //Distance at which bad channels should be serached
949 Int_t relid[4]={0,0,0,0} ;
950 fPHOSGeo->RelPosToRelId(mod, locPos->X(), locPos->Z(), relid) ;
951 Int_t xmin=TMath::Max(1,relid[2]-range) ;
952 Int_t xmax=TMath::Min(64,relid[2]+range) ;
953 Int_t zmin=TMath::Max(1,relid[3]-range) ;
954 Int_t zmax=TMath::Min(56,relid[3]+range) ;
957 for(Int_t ix=xmin;ix<=xmax;ix++){
958 for(Int_t iz=zmin;iz<=zmax;iz++){
959 if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0){ //Bad channel
960 Int_t relidBC[4]={mod,0,ix,iz} ;
961 fPHOSGeo->RelPosInModule(relidBC,x,z);
962 Double_t dist = TMath::Sqrt((x-locPos->X())*(x-locPos->X()) + (z-locPos->Z())*(z-locPos->Z()));
963 if(dist<minDist) minDist = dist;