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 //Correct Misalignment
313 // CorrectPHOSMisalignment(global,mod) ;
314 // position[0]=global.X() ;
315 // position[1]=global.Y() ;
316 // position[2]=global.Z() ;
317 // cluPHOS.GetPosition(position) ;
319 //Eval CoreDispersion
320 Double_t m02=0.,m20=0.;
321 EvalLambdas(&cluPHOS,m02, m20);
322 clu->SetDispersion(TestLambda(clu->E(),m20,m02)) ;
325 cluPHOS.GetPosition(xyz);
326 clu->SetPosition(xyz); //rec.point position in MARS
327 clu->SetE(cluPHOS.E()); //total or core particle energy
328 // clu->SetDispersion(cluPHOS.GetDispersion()); //cluster dispersion
329 // ec->SetPID(rp->GetPID()) ; //array of particle identification
330 clu->SetM02(cluPHOS.GetM02()) ; //second moment M2x
331 clu->SetM20(cluPHOS.GetM20()) ; //second moment M2z
332 Double_t dx=0.,dz=0. ;
333 fPHOSGeo->GlobalPos2RelId(global,relId) ;
335 fPHOSGeo->Global2Local(locPos,global,mod) ;
339 FindTrackMatching(mod,&locPos,dx,dz,pttrack,charge) ;
340 Double_t r=TestCPV(dx, dz, pttrack,charge) ;
341 clu->SetTrackDistance(dx,dz);
343 clu->SetEmcCpvDistance(r);
344 clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02())); //not yet implemented
345 clu->SetTOF(EvalTOF(&cluPHOS,cells));
346 Double_t minDist=clu->GetDistanceToBadChannel() ;//Already calculated
347 DistanceToBadChannel(mod,&locPos,minDist);
348 clu->SetDistanceToBadChannel(minDist) ;
353 Int_t multClust=aod->GetNumberOfCaloClusters();
354 AliAODCaloCells * cells = aod->GetPHOSCells() ;
356 for (Int_t i=0; i<multClust; i++) {
357 AliAODCaloCluster *clu = aod->GetCaloCluster(i);
358 if ( !clu->IsPHOS()) continue;
361 clu->GetPosition(position);
362 TVector3 global(position) ;
364 fPHOSGeo->GlobalPos2RelId(global,relId) ;
365 Int_t mod = relId[0] ;
366 Int_t cellX = relId[2];
367 Int_t cellZ = relId[3] ;
368 if ( !IsGoodChannel(mod,cellX,cellZ) ) {
372 TVector3 locPosOld; //Use it to re-calculate distance to track
373 fPHOSGeo->Global2Local(locPosOld,global,mod) ;
375 //Apply re-Calibreation
376 AliPHOSAodCluster cluPHOS(*clu);
377 cluPHOS.Recalibrate(fPHOSCalibData,cells); // modify the cell energies
378 cluPHOS.EvalAll(logWeight,vertex); // recalculate the cluster parameters
379 cluPHOS.SetE(CorrectNonlinearity(cluPHOS.E()));// Users's nonlinearity
381 //Correct Misalignment
382 // cluPHOS.GetPosition(position);
383 // global.SetXYZ(position[0],position[1],position[2]);
384 // CorrectPHOSMisalignment(global,mod) ;
385 // position[0]=global.X() ;
386 // position[1]=global.Y() ;
387 // position[2]=global.Z() ;
389 clu->SetPosition(position); //rec.point position in MARS
390 clu->SetE(cluPHOS.E()); //total or core particle energy
391 clu->SetDispersion(cluPHOS.GetDispersion()); //cluster dispersion
392 // ec->SetPID(rp->GetPID()) ; //array of particle identification
393 clu->SetM02(cluPHOS.GetM02()) ; //second moment M2x
394 clu->SetM20(cluPHOS.GetM20()) ; //second moment M2z
395 //correct distance to track
396 Double_t dx=clu->GetTrackDx() ;
397 Double_t dz=clu->GetTrackDz() ;
399 fPHOSGeo->Global2Local(locPos,global,mod) ;
400 if(dx!=-999.){ //there is matched track
401 dx+=locPos.X()-locPosOld.X() ;
402 dz+=locPos.Z()-locPosOld.Z() ;
403 clu->SetTrackDistance(dx,dz);
405 Double_t r = 999. ; //Big distance
406 int nTracksMatched = clu->GetNTracksMatched();
407 if(nTracksMatched > 0) {
408 AliVTrack* track = dynamic_cast<AliVTrack*> (clu->GetTrackMatched(0));
410 Double_t pttrack = track->Pt();
411 Short_t charge = track->Charge();
412 r=TestCPV(dx, dz, pttrack,charge) ;
415 clu->SetEmcCpvDistance(r); //Distance in sigmas
417 clu->SetChi2(TestLambda(clu->E(),clu->GetM20(),clu->GetM02())); //not yet implemented
418 clu->SetTOF(EvalTOF(&cluPHOS,cells));
419 Double_t minDist=clu->GetDistanceToBadChannel() ;//Already calculated
420 DistanceToBadChannel(mod,&locPos,minDist);
421 clu->SetDistanceToBadChannel(minDist) ;
426 //___________________________________________________________________________________________________
427 void AliPHOSTenderSupply::FindTrackMatching(Int_t mod,TVector3 *locpos,
428 Double_t &dx, Double_t &dz,
429 Double_t &pt,Int_t &charge){
430 //Find track with closest extrapolation to cluster
431 AliESDEvent *esd = 0x0 ;
433 esd= fTender->GetEvent();
435 esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
439 AliError("ESD is not found") ;
442 Double_t magF = esd->GetMagneticField();
444 Double_t magSign = 1.0;
445 if(magF<0)magSign = -1.0;
447 if (!TGeoGlobalMagField::Instance()->GetField()) {
448 AliMagF* field = new AliMagF("Maps","Maps", magSign, magSign, AliMagF::k5kG);
449 TGeoGlobalMagField::Instance()->SetField(field);
452 // *** Start the matching
454 nt = esd->GetNumberOfTracks();
455 //Calculate actual distance to PHOS module
457 fPHOSGeo->Local2Global(mod, 0.,0., globaPos) ;
458 const Double_t rPHOS = globaPos.Pt() ; //Distance to center of PHOS module
459 const Double_t kYmax = 72.+10. ; //Size of the module (with some reserve) in phi direction
460 const Double_t kZmax = 64.+10. ; //Size of the module (with some reserve) in z direction
461 const Double_t kAlpha0=330./180.*TMath::Pi() ; //First PHOS module angular direction
462 const Double_t kAlpha= 20./180.*TMath::Pi() ; //PHOS module angular size
463 Double_t minDistance = 1.e6;
466 Double_t gposTrack[3] ;
468 Double_t bz = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->SolenoidField();
469 bz = TMath::Sign(0.5*kAlmost0Field,bz) + bz;
473 for (Int_t i=0; i<nt; i++) {
474 AliESDtrack *esdTrack=esd->GetTrack(i);
476 // Skip the tracks having "wrong" status (has to be checked/tuned)
477 ULong_t status = esdTrack->GetStatus();
478 if ((status & AliESDtrack::kTPCout) == 0) continue;
480 //Continue extrapolation from TPC outer surface
481 const AliExternalTrackParam *outerParam=esdTrack->GetOuterParam();
482 if (!outerParam) continue;
483 AliExternalTrackParam t(*outerParam);
486 //Direction to the current PHOS module
487 Double_t phiMod=kAlpha0-kAlpha*mod ;
488 if(!t.Rotate(phiMod))
491 Double_t y; // Some tracks do not reach the PHOS
492 if (!t.GetYAt(rPHOS,bz,y)) continue; // because of the bending
495 if(!t.GetZAt(rPHOS,bz,z))
497 if (TMath::Abs(z) > kZmax)
498 continue; // Some tracks miss the PHOS in Z
499 if(TMath::Abs(y) < kYmax){
500 t.PropagateToBxByBz(rPHOS,b); // Propagate to the matching module
501 //t.CorrectForMaterial(...); // Correct for the TOF material, if needed
502 t.GetXYZ(gposTrack) ;
503 TVector3 globalPositionTr(gposTrack) ;
504 TVector3 localPositionTr ;
505 fPHOSGeo->Global2Local(localPositionTr,globalPositionTr,mod) ;
506 Double_t ddx = locpos->X()-localPositionTr.X();
507 Double_t ddz = locpos->Z()-localPositionTr.Z();
508 Double_t d2 = ddx*ddx + ddz*ddz;
509 if(d2 < minDistance) {
514 charge=esdTrack->Charge() ;
517 }//Scanned all tracks
520 //____________________________________________________________
521 Float_t AliPHOSTenderSupply::CorrectNonlinearity(Float_t en){
523 //For backward compatibility, if no RecoParameters found
524 if(fNonlinearityVersion=="Default"){
525 return 0.0241+1.0504*en+0.000249*en*en ;
528 if(fNonlinearityVersion=="NoCorrection"){
531 if(fNonlinearityVersion=="Gustavo2005"){
532 return fNonlinearityParams[0]+fNonlinearityParams[1]*en + fNonlinearityParams[2]*en*en ;
534 if(fNonlinearityVersion=="Henrik2010"){
535 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])) ;
540 //_____________________________________________________________________________
541 Double_t AliPHOSTenderSupply::TestLambda(Double_t pt,Double_t l1,Double_t l2){
542 //Parameterization for core dispersion
544 Double_t l1Mean = 1.150200 + 0.097886/(1.+1.486645*pt+0.000038*pt*pt) ;
545 Double_t l2Mean = 1.574706 + 0.997966*exp(-0.895075*pt)-0.010666*pt ;
546 Double_t l1Sigma = 0.100255 + 0.337177*exp(-0.517684*pt)+0.001170*pt ;
547 Double_t l2Sigma = 0.232580 + 0.573401*exp(-0.735903*pt)-0.002325*pt ;
548 Double_t c = -0.110983 -0.017353/(1.-1.836995*pt+0.934517*pt*pt) ;
551 //Parameterizatino for full dispersion
552 Double_t l2Mean = 1.53126+9.50835e+06/(1.+1.08728e+07*pt+1.73420e+06*pt*pt) ;
553 Double_t l1Mean = 1.12365+0.123770*TMath::Exp(-pt*0.246551)+5.30000e-03*pt ;
554 Double_t l2Sigma = 6.48260e-02+7.60261e+10/(1.+1.53012e+11*pt+5.01265e+05*pt*pt)+9.00000e-03*pt;
555 Double_t l1Sigma = 4.44719e-04+6.99839e-01/(1.+1.22497e+00*pt+6.78604e-07*pt*pt)+9.00000e-03*pt;
556 Double_t c=-0.35-0.550*TMath::Exp(-0.390730*pt) ;
558 Double_t R2=0.5*(l1-l1Mean)*(l1-l1Mean)/l1Sigma/l1Sigma +
559 0.5*(l2-l2Mean)*(l2-l2Mean)/l2Sigma/l2Sigma +
560 0.5*c*(l1-l1Mean)*(l2-l2Mean)/l1Sigma/l2Sigma ;
564 //____________________________________________________________________________
565 Double_t AliPHOSTenderSupply::TestCPV(Double_t dx, Double_t dz, Double_t pt, Int_t charge){
566 //Parameterization of LHC10h period
571 Double_t sx=TMath::Min(5.4,2.59719e+02*TMath::Exp(-pt/1.02053e-01)+
572 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);
573 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) ;
575 Double_t mf = 0.; //Positive for ++ and negative for --
577 AliESDEvent *esd = fTender->GetEvent();
578 mf = esd->GetMagneticField();
582 AliESDEvent *esd= dynamic_cast<AliESDEvent*>(fTask->InputEvent());
584 mf = esd->GetMagneticField();
586 AliAODEvent *aod= dynamic_cast<AliAODEvent*>(fTask->InputEvent());
588 mf = aod->GetMagneticField();
591 AliError("Neither Tender nor Task defined") ;
595 if(mf<0.){ //field --
598 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)) ;
600 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)) ;
605 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)) ;
607 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)) ;
610 Double_t rz=(dz-meanZ)/sz ;
611 Double_t rx=(dx-meanX)/sx ;
612 return TMath::Sqrt(rx*rx+rz*rz) ;
615 //________________________________________________________________________
616 Bool_t AliPHOSTenderSupply::IsGoodChannel(Int_t mod, Int_t ix, Int_t iz)
618 //Check if this channel belogs to the good ones
621 // AliError(Form("No bad map for PHOS module %d ",mod)) ;
624 if(!fPHOSBadMap[mod]){
625 // AliError(Form("No Bad map for PHOS module %d",mod)) ;
628 if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0)
633 //________________________________________________________________________
634 void AliPHOSTenderSupply::ForceUsingBadMap(const char * filename){
635 //Read TH2I histograms with bad maps from local or alien file
636 TFile * fbm = TFile::Open(filename) ;
637 if(!fbm || !fbm->IsOpen()){
638 AliError(Form("Can not open BadMaps file %s",filename)) ;
643 for(Int_t mod=1;mod<4; mod++){
644 snprintf(key,55,"PHOS_BadMap_mod%d",mod) ;
645 TH2I * h = (TH2I*)fbm->Get(key) ;
647 fPHOSBadMap[mod] = new TH2I(*h) ;
650 fUsePrivateBadMap=kTRUE ;
652 //________________________________________________________________________
653 void AliPHOSTenderSupply::ForceUsingCalibration(const char * filename){
654 //Read PHOS recalibration parameters from the file.
655 //We assume that file contains single entry: AliPHOSCalibData
656 TFile * fc = TFile::Open(filename) ;
657 if(!fc || !fc->IsOpen()){
658 AliFatal(Form("Can not open Calibration file %s",filename)) ;
661 fPHOSCalibData = (AliPHOSCalibData*)fc->Get("PHOSCalibration") ;
663 fUsePrivateCalib=kTRUE;
665 //________________________________________________________________________
666 void AliPHOSTenderSupply::CorrectPHOSMisalignment(TVector3 &global,Int_t mod){
667 //Correct for PHOS modules misalignment
669 //correct misalignment
670 const Float_t shiftX[6]={0.,-2.3,-2.11,-1.53,0.,0.} ;
671 const Float_t shiftZ[6]={0.,-0.4, 0.52, 0.8,0.,0.} ;
673 fPHOSGeo->Global2Local(localPos,global,mod) ;
674 fPHOSGeo->Local2Global(mod,localPos.X()+shiftX[mod],localPos.Z()+shiftZ[mod],global);
676 //________________________________________________________________________
677 void AliPHOSTenderSupply::EvalLambdas(AliVCluster * clu, Double_t &m02, Double_t &m20){
678 //calculate dispecrsion of the cluster in the circle with radius distanceCut around the maximum
680 const Double_t rCut=4.5 ;
682 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
683 // Calculates the center of gravity in the local PHOS-module coordinates
685 Double_t xc[100]={0} ;
686 Double_t zc[100]={0} ;
689 Int_t mulDigit=TMath::Min(100,clu->GetNCells()) ;
690 const Double_t logWeight=4.5 ;
691 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
695 fPHOSGeo->AbsToRelNumbering(clu->GetCellAbsId(iDigit), relid) ;
696 fPHOSGeo->RelPosInModule(relid, xi, zi);
699 if (clu->E()>0 && elist[iDigit]>0) {
700 Float_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
701 x += xc[iDigit] * w ;
702 z += zc[iDigit] * w ;
717 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
718 if (clu->E()>0 && elist[iDigit]>0.) {
719 Double_t w = TMath::Max( 0., logWeight + TMath::Log( elist[iDigit] / clu->E() ) ) ;
720 Double_t xi= xc[iDigit] ;
721 Double_t zi= zc[iDigit] ;
722 if((xi-x)*(xi-x)+(zi-z)*(zi-z) < rCut*rCut){
743 m02 = 0.5 * (dxx + dzz) + TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
744 m20 = 0.5 * (dxx + dzz) - TMath::Sqrt( 0.25 * (dxx - dzz) * (dxx - dzz) + dxz * dxz ) ;
751 //________________________________________________________________________
752 Double_t AliPHOSTenderSupply::EvalTOF(AliVCluster * clu,AliVCaloCells * cells){
753 //Evaluate TOF of the cluster after re-calibration
754 //TOF here is weighted average of digits
755 // -within 50ns from the most energetic cell
759 Double32_t * elist = clu->GetCellsAmplitudeFraction() ;
760 Int_t mulDigit=clu->GetNCells() ;
762 Float_t tMax= 0.; //Time at the maximum
764 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
765 Int_t absId=clu->GetCellAbsId(iDigit) ;
767 if(cells->GetCellMCLabel(absId)==-2) //This is LG digit. No statistics to calibrate LG timing, remove them from TOF calculation
769 if( elist[iDigit]>eMax){
770 tMax=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
776 //Try to improve accuracy
777 //Do not account time of soft cells:
778 // const Double_t part=0.5 ;
779 Double_t eMin=TMath::Min(0.5,0.2*eMax) ;
782 for(Int_t iDigit=0; iDigit<mulDigit; iDigit++) {
783 Int_t absId=clu->GetCellAbsId(iDigit) ;
785 if(cells->GetCellMCLabel(absId)==-2) //This is LG digit. No statistics to calibrate LG timing, remove them from TOF calculation
788 Double_t ti=CalibrateTOF(cells->GetCellTime(absId),absId,isHG) ;
789 if(TMath::Abs(ti-tMax)>50.e-9) //remove soft cells with wrong time
792 //Remove too soft cells
793 if(elist[iDigit]<eMin)
797 //weight = 1./sigma^2
798 //Sigma is parameterization of TOF resolution 16.05.2013
801 wi2=1./(2.4e-9 + 3.9e-9/elist[iDigit]) ;
803 wi2=1./(2.4e-9 + 3.9e-9/(0.1*elist[iDigit])) ; //E of LG digit is 1/16 of correcponding HG
815 //________________________________________________________________________
816 Double_t AliPHOSTenderSupply::CalibrateTOF(Double_t tof, Int_t absId, Bool_t isHG){
817 //Apply time re-calibration separately for HG and LG channels
818 //By default (if not filled) shifts are zero.
821 fPHOSGeo->AbsToRelNumbering(absId,relId) ;
822 Int_t module = relId[0];
823 Int_t column = relId[3];
824 Int_t row = relId[2];
826 tof-=fPHOSCalibData->GetTimeShiftEmc(module, column, row);
828 tof-=fPHOSCalibData->GetLGTimeShiftEmc(module, column, row);
833 //________________________________________________________________________
834 void AliPHOSTenderSupply::DistanceToBadChannel(Int_t mod, TVector3 * locPos, Double_t &minDist){
835 //Check if distance to bad channel was reduced
836 Int_t range = minDist/2.2 +1 ; //Distance at which bad channels should be serached
838 Int_t relid[4]={0,0,0,0} ;
839 fPHOSGeo->RelPosToRelId(mod, locPos->X(), locPos->Z(), relid) ;
840 Int_t xmin=TMath::Max(1,relid[2]-range) ;
841 Int_t xmax=TMath::Min(64,relid[2]+range) ;
842 Int_t zmin=TMath::Max(1,relid[3]-range) ;
843 Int_t zmax=TMath::Min(56,relid[3]+range) ;
846 for(Int_t ix=xmin;ix<=xmax;ix++){
847 for(Int_t iz=zmin;iz<=zmax;iz++){
848 if(fPHOSBadMap[mod]->GetBinContent(ix,iz)>0){ //Bad channel
849 Int_t relidBC[4]={mod,0,ix,iz} ;
850 fPHOSGeo->RelPosInModule(relidBC,x,z);
851 Double_t dist = TMath::Sqrt((x-locPos->X())*(x-locPos->X()) + (z-locPos->Z())*(z-locPos->Z()));
852 if(dist<minDist) minDist = dist;