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 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() ;
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");
179 AliError(Form("Can not read Bad map for run %d. \n You may choose to use your map with ForceUsingBadMap()\n",runNumber)) ;
182 AliInfo(Form("Setting PHOS bad map with name %s \n",maps->GetName())) ;
183 for(Int_t mod=0; mod<5;mod++){
185 delete fPHOSBadMap[mod] ;
186 TH2I * h = (TH2I*)maps->At(mod) ;
188 fPHOSBadMap[mod]=new TH2I(*h) ;
193 if(!fUsePrivateCalib){
194 if(fIsMC){ //re/de-calibration for MC productions
196 AliOADBContainer calibContainer("phosRecalibration");
197 calibContainer.InitFromFile("$ALICE_ROOT/OADB/PHOS/PHOSMCCalibrations.root","phosRecalibration");
199 TObjArray *recalib = (TObjArray*)calibContainer.GetObject(runNumber,"PHOSRecalibration");
201 AliFatal(Form("Can not read calibrations for run %d\n. You may choose your specific calibration with ForceUsingCalibration()\n",runNumber)) ;
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 ;
212 if(!fPHOSCalibData) {
213 AliFatal(Form("Can not find calibration for run %d, and name %s \n",runNumber, fMCProduction.Data())) ;
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");
225 AliFatal(Form("Can not read calibrations for run %d\n. You may choose your specific calibration with ForceUsingCalibration()\n",runNumber)) ;
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)) ;
238 //_____________________________________________________
239 void AliPHOSTenderSupply::ProcessEvent()
241 //Choose PHOS clusters and recalibrate them
242 //that it recalculate energy, position and distance
243 //to closest track extrapolation
245 AliESDEvent *esd = 0x0 ;
246 AliAODEvent *aod = 0x0 ;
248 esd = fTender->GetEvent();
256 esd = dynamic_cast<AliESDEvent*>(fTask->InputEvent()) ;
257 aod = dynamic_cast<AliAODEvent*>(fTask->InputEvent()) ;
263 || (fTender && fTender->RunChanged())){ //In case of Task init called automatically
270 const AliESDVertex *esdVertex = esd->GetPrimaryVertex();
271 vertex.SetXYZ(esdVertex->GetX(),esdVertex->GetY(),esdVertex->GetZ());
274 const AliAODVertex *aodVertex = aod->GetPrimaryVertex();
275 vertex.SetXYZ(aodVertex->GetX(),aodVertex->GetY(),aodVertex->GetZ());
277 if(vertex.Mag()>99.) //vertex not defined?
278 vertex.SetXYZ(0.,0.,0.) ;
282 const Double_t logWeight=4.5 ;
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() ;
289 for (Int_t i=0; i<multClust; i++) {
290 AliESDCaloCluster *clu = esd->GetCaloCluster(i);
291 if ( !clu->IsPHOS()) continue;
294 clu->GetPosition(position);
295 TVector3 global(position) ;
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) ) {
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
312 Double_t ecore=CoreEnergy(&cluPHOS) ;
313 ecore=CorrectNonlinearity(ecore) ;
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) ;
322 //Eval CoreDispersion
323 Double_t m02=0.,m20=0.;
324 EvalLambdas(&cluPHOS,m02, m20);
325 clu->SetDispersion(TestLambda(clu->E(),m20,m02)) ;
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
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) ;
340 fPHOSGeo->Global2Local(locPos,global,mod) ;
344 FindTrackMatching(mod,&locPos,dx,dz,pttrack,charge) ;
345 Double_t r=TestCPV(dx, dz, pttrack,charge) ;
346 clu->SetTrackDistance(dx,dz);
348 clu->SetEmcCpvDistance(r);
349 clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02())); //not yet implemented
350 Double_t tof=EvalTOF(&cluPHOS,cells);
351 // if(TMath::Abs(tof-clu->GetTOF())>100.e-9) //something wrong in cell TOF!
352 // tof=clu->GetTOF() ;
354 Double_t minDist=clu->GetDistanceToBadChannel() ;//Already calculated
355 DistanceToBadChannel(mod,&locPos,minDist);
356 clu->SetDistanceToBadChannel(minDist) ;
361 Int_t multClust=aod->GetNumberOfCaloClusters();
362 AliAODCaloCells * cells = aod->GetPHOSCells() ;
364 for (Int_t i=0; i<multClust; i++) {
365 AliAODCaloCluster *clu = aod->GetCaloCluster(i);
366 if ( !clu->IsPHOS()) continue;
369 clu->GetPosition(position);
370 TVector3 global(position) ;
372 fPHOSGeo->GlobalPos2RelId(global,relId) ;
373 Int_t mod = relId[0] ;
374 Int_t cellX = relId[2];
375 Int_t cellZ = relId[3] ;
376 if ( !IsGoodChannel(mod,cellX,cellZ) ) {
380 TVector3 locPosOld; //Use it to re-calculate distance to track
381 fPHOSGeo->Global2Local(locPosOld,global,mod) ;
383 //Apply re-Calibreation
384 AliPHOSAodCluster cluPHOS(*clu);
385 cluPHOS.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
386 cluPHOS.EvalAll(logWeight,vertex); // recalculate the cluster parameters
387 cluPHOS.SetE(CorrectNonlinearity(cluPHOS.E()));// Users's nonlinearity
389 Double_t ecore=CoreEnergy(&cluPHOS) ;
390 ecore=CorrectNonlinearity(ecore) ;
393 //Correct Misalignment
394 // cluPHOS.GetPosition(position);
395 // global.SetXYZ(position[0],position[1],position[2]);
396 // CorrectPHOSMisalignment(global,mod) ;
397 // position[0]=global.X() ;
398 // position[1]=global.Y() ;
399 // position[2]=global.Z() ;
401 clu->SetPosition(position); //rec.point position in MARS
402 clu->SetE(cluPHOS.E()); //total particle energy
403 clu->SetCoreEnergy(ecore); //core particle energy
404 clu->SetDispersion(cluPHOS.GetDispersion()); //cluster dispersion
405 // ec->SetPID(rp->GetPID()) ; //array of particle identification
406 clu->SetM02(cluPHOS.GetM02()) ; //second moment M2x
407 clu->SetM20(cluPHOS.GetM20()) ; //second moment M2z
408 //correct distance to track
409 Double_t dx=clu->GetTrackDx() ;
410 Double_t dz=clu->GetTrackDz() ;
412 fPHOSGeo->Global2Local(locPos,global,mod) ;
413 if(dx!=-999.){ //there is matched track
414 dx+=locPos.X()-locPosOld.X() ;
415 dz+=locPos.Z()-locPosOld.Z() ;
416 clu->SetTrackDistance(dx,dz);
418 Double_t r = 999. ; //Big distance
419 int nTracksMatched = clu->GetNTracksMatched();
420 if(nTracksMatched > 0) {
421 AliVTrack* track = dynamic_cast<AliVTrack*> (clu->GetTrackMatched(0));
423 Double_t pttrack = track->Pt();
424 Short_t charge = track->Charge();
425 r=TestCPV(dx, dz, pttrack,charge) ;
428 clu->SetEmcCpvDistance(r); //Distance in sigmas
430 clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02())); //not yet implemented
431 Double_t tof=EvalTOF(&cluPHOS,cells);
432 // if(TMath::Abs(tof-clu->GetTOF())>100.e-9) //something wrong in cell TOF!
433 // tof=clu->GetTOF() ;
435 Double_t minDist=clu->GetDistanceToBadChannel() ;//Already calculated
436 DistanceToBadChannel(mod,&locPos,minDist);
437 clu->SetDistanceToBadChannel(minDist) ;
438 Double_t eCross=EvalEcross(&cluPHOS);
439 clu->SetMCEnergyFraction(eCross) ;
444 //___________________________________________________________________________________________________
445 void AliPHOSTenderSupply::FindTrackMatching(Int_t mod,TVector3 *locpos,
446 Double_t &dx, Double_t &dz,
447 Double_t &pt,Int_t &charge){
448 //Find track with closest extrapolation to cluster
449 AliESDEvent *esd = 0x0 ;
451 esd= fTender->GetEvent();
453 esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
457 AliError("ESD is not found") ;
460 Double_t magF = esd->GetMagneticField();
462 Double_t magSign = 1.0;
463 if(magF<0)magSign = -1.0;
465 if (!TGeoGlobalMagField::Instance()->GetField()) {
466 AliError("Margnetic filed was not initialized, use default") ;
467 AliMagF* field = new AliMagF("Maps","Maps", magSign, magSign, AliMagF::k5kG);
468 TGeoGlobalMagField::Instance()->SetField(field);
471 // *** Start the matching
473 nt = esd->GetNumberOfTracks();
474 //Calculate actual distance to PHOS module
476 fPHOSGeo->Local2Global(mod, 0.,0., globaPos) ;
477 const Double_t rPHOS = globaPos.Pt() ; //Distance to center of PHOS module
478 const Double_t kYmax = 72.+10. ; //Size of the module (with some reserve) in phi direction
479 const Double_t kZmax = 64.+10. ; //Size of the module (with some reserve) in z direction
480 const Double_t kAlpha0=330./180.*TMath::Pi() ; //First PHOS module angular direction
481 const Double_t kAlpha= 20./180.*TMath::Pi() ; //PHOS module angular size
482 Double_t minDistance = 1.e6;
485 Double_t gposTrack[3] ;
487 Double_t bz = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField();
488 bz = TMath::Sign(0.5*kAlmost0Field,bz) + bz;
492 for (Int_t i=0; i<nt; i++) {
493 AliESDtrack *esdTrack=esd->GetTrack(i);
495 // Skip the tracks having "wrong" status (has to be checked/tuned)
496 ULong_t status = esdTrack->GetStatus();
497 if ((status & AliESDtrack::kTPCout) == 0) continue;
499 //Continue extrapolation from TPC outer surface
500 const AliExternalTrackParam *outerParam=esdTrack->GetOuterParam();
501 if (!outerParam) continue;
502 AliExternalTrackParam t(*outerParam);
505 //Direction to the current PHOS module
506 Double_t phiMod=kAlpha0-kAlpha*mod ;
507 if(!t.Rotate(phiMod))
510 Double_t y; // Some tracks do not reach the PHOS
511 if (!t.GetYAt(rPHOS,bz,y)) continue; // because of the bending
514 if(!t.GetZAt(rPHOS,bz,z))
516 if (TMath::Abs(z) > kZmax)
517 continue; // Some tracks miss the PHOS in Z
518 if(TMath::Abs(y) < kYmax){
519 t.PropagateToBxByBz(rPHOS,b); // Propagate to the matching module
520 //t.CorrectForMaterial(...); // Correct for the TOF material, if needed
521 t.GetXYZ(gposTrack) ;
522 TVector3 globalPositionTr(gposTrack) ;
523 TVector3 localPositionTr ;
524 fPHOSGeo->Global2Local(localPositionTr,globalPositionTr,mod) ;
525 Double_t ddx = locpos->X()-localPositionTr.X();
526 Double_t ddz = locpos->Z()-localPositionTr.Z();
527 Double_t d2 = ddx*ddx + ddz*ddz;
528 if(d2 < minDistance) {
533 charge=esdTrack->Charge() ;
536 }//Scanned all tracks
539 //____________________________________________________________
540 Float_t AliPHOSTenderSupply::CorrectNonlinearity(Float_t en){
542 //For backward compatibility, if no RecoParameters found
543 if(fNonlinearityVersion=="Default"){
544 return 0.0241+1.0504*en+0.000249*en*en ;
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::TestLambda(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::TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge){
585 //Parameterization of LHC10h period
590 Double_t sx=TMath::Min(5.4,2.59719e+02*TMath::Exp(-pt/1.02053e-01)+
591 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);
592 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) ;
594 Double_t mf = 0.; //Positive for ++ and negative for --
596 AliESDEvent *esd = fTender->GetEvent();
597 mf = esd->GetMagneticField();
601 AliESDEvent *esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
603 mf = esd->GetMagneticField();
605 AliAODEvent *aod= dynamic_cast<AliAODEvent*>(fTask->InputEvent());
607 mf = aod->GetMagneticField();
610 AliError("Neither Tender nor Task defined") ;
614 if(mf<0.){ //field --
617 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)) ;
619 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)) ;
624 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)) ;
626 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)) ;
629 Double_t rz=(dz-meanZ)/sz ;
630 Double_t rx=(dx-meanX)/sx ;
631 return TMath::Sqrt(rx*rx+rz*rz) ;
634 //________________________________________________________________________
635 Bool_t AliPHOSTenderSupply::IsGoodChannel(Int_t mod, Int_t ix, Int_t iz)
637 //Check if this channel belogs to the good ones
640 // AliError(Form("No bad map for PHOS module %d ",mod)) ;
643 if(!fPHOSBadMap[mod]){
644 // AliError(Form("No Bad map for PHOS module %d",mod)) ;
647 if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0)
652 //________________________________________________________________________
653 void AliPHOSTenderSupply::ForceUsingBadMap(const char * filename){
654 //Read TH2I histograms with bad maps from local or alien file
655 TFile * fbm = TFile::Open(filename) ;
656 if(!fbm || !fbm->IsOpen()){
657 AliError(Form("Can not open BadMaps file %s",filename)) ;
662 for(Int_t mod=1;mod<4; mod++){
663 snprintf(key,55,"PHOS_BadMap_mod%d",mod) ;
664 TH2I * h = (TH2I*)fbm->Get(key) ;
666 fPHOSBadMap[mod] = new TH2I(*h) ;
669 fUsePrivateBadMap=kTRUE ;
671 //________________________________________________________________________
672 void AliPHOSTenderSupply::ForceUsingCalibration(const char * filename){
673 //Read PHOS recalibration parameters from the file.
674 //We assume that file contains single entry: AliPHOSCalibData
675 TFile * fc = TFile::Open(filename) ;
676 if(!fc || !fc->IsOpen()){
677 AliFatal(Form("Can not open Calibration file %s",filename)) ;
680 fPHOSCalibData = (AliPHOSCalibData*)fc->Get("PHOSCalibration") ;
682 fUsePrivateCalib=kTRUE;
684 //________________________________________________________________________
685 void AliPHOSTenderSupply::CorrectPHOSMisalignment(TVector3 &global,Int_t mod){
686 //Correct for PHOS modules misalignment
688 //correct misalignment
689 const Float_t shiftX[6]={0.,-2.3,-2.11,-1.53,0.,0.} ;
690 const Float_t shiftZ[6]={0.,-0.4, 0.52, 0.8,0.,0.} ;
692 fPHOSGeo->Global2Local(localPos,global,mod) ;
693 fPHOSGeo->Local2Global(mod,localPos.X()+shiftX[mod],localPos.Z()+shiftZ[mod],global);
695 //________________________________________________________________________
696 void AliPHOSTenderSupply::EvalLambdas(AliVCluster * clu, Double_t &m02, Double_t &m20){
697 //calculate dispecrsion of the cluster in the circle with radius distanceCut around the maximum
699 const Double_t rCut=4.5 ;
701 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
702 // Calculates the center of gravity in the local PHOS-module coordinates
704 Double_t xc[100]={0} ;
705 Double_t zc[100]={0} ;
708 Int_t mulDigit=TMath::Min(100,clu->GetNCells()) ;
709 const Double_t logWeight=4.5 ;
710 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
714 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
715 fPHOSGeo->RelPosInModule(relid, xi, zi);
718 if (clu->E()>0 && elist[iDigit]>0) {
719 Float_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
720 x += xc[iDigit] * w ;
721 z += zc[iDigit] * w ;
736 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
737 if (clu->E()>0 && elist[iDigit]>0.) {
738 Double_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
739 Double_t xi= xc[iDigit] ;
740 Double_t zi= zc[iDigit] ;
741 if((xi-x)*(xi-x)+(zi-z)*(zi-z) < rCut*rCut){
762 m02 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
763 m20 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
770 //____________________________________________________________________________
771 Double_t AliPHOSTenderSupply::CoreEnergy(AliVCluster * clu){
772 //calculate energy of the cluster in the circle with radius distanceCut around the maximum
774 //Can not use already calculated coordinates?
775 //They have incidence correction...
776 const Double_t distanceCut =3.5 ;
777 const Double_t logWeight=4.5 ;
779 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
780 // Calculates the center of gravity in the local PHOS-module coordinates
782 Double_t xc[100]={0} ;
783 Double_t zc[100]={0} ;
786 Int_t mulDigit=TMath::Min(100,clu->GetNCells()) ;
787 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
791 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
792 fPHOSGeo->RelPosInModule(relid, xi, zi);
795 if (clu->E()>0 && elist[iDigit]>0) {
796 Float_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
797 x += xc[iDigit] * w ;
798 z += zc[iDigit] * w ;
807 for(Int_t iDigit=0; iDigit < mulDigit; iDigit++) {
808 Double_t distance = TMath::Sqrt((xc[iDigit]-x)*(xc[iDigit]-x)+(zc[iDigit]-z)*(zc[iDigit]-z)) ;
809 if(distance < distanceCut)
810 coreE += elist[iDigit] ;
812 //Apply non-linearity correction
815 //____________________________________________________________________________
816 Double_t AliPHOSTenderSupply::EvalEcross(AliVCluster * clu){
817 //Calculate propoerion of the cluster energy in cross around the
818 //cell with maximal energy deposition. Can be used to reject exotic clusters
820 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
821 Int_t mulDigit=clu->GetNCells() ;
822 // Calculates the center of gravity in the local PHOS-module coordinates
823 //Find cell with max E
826 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
827 if(elist[iDigit]>eMax){
832 //Calculate e in cross
835 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iMax), relidMax) ;
837 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
838 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
839 if(TMath::Abs(relid[2]-relidMax[2])+TMath::Abs(relid[3]-relidMax[3])==1)
840 eCross+= elist[iDigit] ;
843 return 1.-eCross/eMax ;
849 //________________________________________________________________________
850 Double_t AliPHOSTenderSupply::EvalTOF(AliVCluster * clu,AliVCaloCells * cells){
851 //Evaluate TOF of the cluster after re-calibration
852 //TOF here is weighted average of digits
853 // -within 50ns from the most energetic cell
856 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
857 Int_t mulDigit=clu->GetNCells() ;
859 Float_t tMax= 0.; //Time at the maximum
861 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
862 Int_t absId=clu->GetCellAbsId(iDigit) ;
864 if(cells->GetCellMCLabel(absId)==-2){ //This is LG digit.
867 if( elist[iDigit]>eMax){
868 tMax=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
873 //Try to improve accuracy
874 //Do not account time of soft cells:
875 // const Double_t part=0.5 ;
876 Double_t eMin=TMath::Min(0.5,0.2*eMax) ;
879 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
880 Int_t absId=clu->GetCellAbsId(iDigit) ;
882 if(cells->GetCellMCLabel(absId)==-2){ //This is LG digit.
886 Double_t ti=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
887 if(TMath::Abs(ti-tMax)>50.e-9) //remove soft cells with wrong time
890 //Remove too soft cells
891 if(elist[iDigit]<eMin)
895 //weight = 1./sigma^2
896 //Sigma is parameterization of TOF resolution 16.05.2013
899 wi2=1./(2.4 + 3.9/elist[iDigit]) ;
901 wi2=1./(2.4 + 3.9/(0.1*elist[iDigit])) ; //E of LG digit is 1/16 of correcponding HG
916 //________________________________________________________________________
917 Double_t AliPHOSTenderSupply::CalibrateTOF(Double_t tof, Int_t absId, Bool_t isHG){
918 //Apply time re-calibration separately for HG and LG channels
919 //By default (if not filled) shifts are zero.
922 fPHOSGeo->AbsToRelNumbering(absId,relId) ;
923 Int_t module = relId[0];
924 Int_t column = relId[3];
925 Int_t row = relId[2];
927 tof-=fPHOSCalibData->GetTimeShiftEmc(module, column, row);
929 tof-=fPHOSCalibData->GetLGTimeShiftEmc(module, column, row);
934 //________________________________________________________________________
935 void AliPHOSTenderSupply::DistanceToBadChannel(Int_t mod, TVector3 * locPos, Double_t &minDist){
936 //Check if distance to bad channel was reduced
937 Int_t range = minDist/2.2 +1 ; //Distance at which bad channels should be serached
939 Int_t relid[4]={0,0,0,0} ;
940 fPHOSGeo->RelPosToRelId(mod, locPos->X(), locPos->Z(), relid) ;
941 Int_t xmin=TMath::Max(1,relid[2]-range) ;
942 Int_t xmax=TMath::Min(64,relid[2]+range) ;
943 Int_t zmin=TMath::Max(1,relid[3]-range) ;
944 Int_t zmax=TMath::Min(56,relid[3]+range) ;
947 for(Int_t ix=xmin;ix<=xmax;ix++){
948 for(Int_t iz=zmin;iz<=zmax;iz++){
949 if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0){ //Bad channel
950 Int_t relidBC[4]={mod,0,ix,iz} ;
951 fPHOSGeo->RelPosInModule(relidBC,x,z);
952 Double_t dist = TMath::Sqrt((x-locPos->X())*(x-locPos->X()) + (z-locPos->Z())*(z-locPos->Z()));
953 if(dist<minDist) minDist = dist;