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 **************************************************************************/
18 //_________________________________________________________________________
20 //-- Yves Schutz (SUBATECH)
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and
22 // derived from STEER/AliReconstructor.
26 // --- ROOT system ---
27 #include <TClonesArray.h>
28 #include "TGeoManager.h"
29 #include "TGeoMatrix.h"
31 // --- Standard library ---
33 // --- AliRoot header files ---
34 #include "AliEMCALReconstructor.h"
36 #include "AliCodeTimer.h"
37 #include "AliCaloCalibPedestal.h"
38 #include "AliEMCALCalibData.h"
39 #include "AliESDEvent.h"
40 #include "AliESDCaloCluster.h"
41 #include "AliESDCaloCells.h"
42 #include "AliESDtrack.h"
43 #include "AliEMCALLoader.h"
44 #include "AliEMCALRawUtils.h"
45 #include "AliEMCALDigit.h"
46 #include "AliEMCALClusterizerv1.h"
47 #include "AliEMCALClusterizerNxN.h"
48 #include "AliEMCALRecPoint.h"
49 #include "AliEMCALPID.h"
50 #include "AliEMCALTrigger.h"
51 #include "AliRawReader.h"
52 #include "AliCDBEntry.h"
53 #include "AliCDBManager.h"
54 #include "AliEMCALGeometry.h"
56 #include "AliESDVZERO.h"
57 #include "AliCDBManager.h"
58 #include "AliRunLoader.h"
60 #include "AliEMCALTriggerData.h"
61 #include "AliEMCALTriggerElectronics.h"
62 #include "AliEMCALTriggerDCSConfigDB.h"
63 #include "AliEMCALTriggerDCSConfig.h"
64 #include "AliEMCALTriggerData.h"
65 #include "AliEMCALTriggerRawDigit.h"
66 #include "AliEMCALTriggerPatch.h"
67 #include "AliEMCALTriggerTypes.h"
69 ClassImp(AliEMCALReconstructor)
71 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0; // EMCAL rec. parameters
72 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0; // EMCAL raw utilities class
73 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0; // EMCAL clusterizer class
74 TClonesArray* AliEMCALReconstructor::fgDigitsArr = 0; // list of digits, to be used multiple times
75 TObjArray* AliEMCALReconstructor::fgClustersArr = 0; // list of clusters, to be used multiple times
76 TClonesArray* AliEMCALReconstructor::fgTriggerDigits = 0; // list of trigger digits, to be used multiple times
77 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
78 //____________________________________________________________________________
79 AliEMCALReconstructor::AliEMCALReconstructor()
80 : fGeom(0),fCalibData(0),fPedestalData(0),fTriggerData(0x0), fMatches(0x0)
84 // AliDebug(2, "Mark.");
86 fgRawUtils = new AliEMCALRawUtils;
88 //To make sure we match with the geometry in a simulation file,
89 //let's try to get it first. If not, take the default geometry
90 AliRunLoader *rl = AliRunLoader::Instance();
92 AliEMCAL * emcal = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"));
93 if(emcal) fGeom = emcal->GetGeometry();
97 AliInfo(Form("Using default geometry in reconstruction"));
98 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
101 //Get calibration parameters
104 AliCDBEntry *entry = (AliCDBEntry*)
105 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
106 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
110 AliFatal("Calibration parameters not found in CDB!");
112 //Get calibration parameters
115 AliCDBEntry *entry = (AliCDBEntry*)
116 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
117 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
121 AliFatal("Dead map not found in CDB!");
123 if(!fGeom) AliFatal(Form("Could not get geometry!"));
125 AliEMCALTriggerDCSConfigDB* dcsConfigDB = AliEMCALTriggerDCSConfigDB::Instance();
127 const AliEMCALTriggerDCSConfig* dcsConfig = dcsConfigDB->GetTriggerDCSConfig();
129 if (!dcsConfig) AliFatal("No Trigger DCS Configuration from OCDB!");
130 fgTriggerProcessor = new AliEMCALTriggerElectronics( dcsConfig );
132 fTriggerData = new AliEMCALTriggerData();
134 //Init temporary list of digits
135 fgDigitsArr = new TClonesArray("AliEMCALDigit",1000);
136 fgClustersArr = new TObjArray(1000);
137 fgTriggerDigits = new TClonesArray("AliEMCALTriggerRawDigit",1000);
140 fMatches = new TList();
141 fMatches->SetOwner(kTRUE);
144 //____________________________________________________________________________
145 AliEMCALReconstructor::~AliEMCALReconstructor()
149 //AliDebug(2, "Mark.");
151 if(fGeom) delete fGeom;
153 //No need to delete, recovered from OCDB
154 //if(fCalibData) delete fCalibData;
155 //if(fPedestalData) delete fPedestalData;
158 fgDigitsArr->Clear("C");
163 fgClustersArr->Clear();
164 delete fgClustersArr;
168 fgTriggerDigits->Clear();
169 delete fgTriggerDigits;
172 if(fgRawUtils) delete fgRawUtils;
173 if(fgClusterizer) delete fgClusterizer;
174 if(fgTriggerProcessor) delete fgTriggerProcessor;
176 if(fMatches) { fMatches->Delete(); delete fMatches; fMatches = 0;}
178 AliCodeTimer::Instance()->Print();
181 //____________________________________________________________________________
182 void AliEMCALReconstructor::InitClusterizer() const
184 //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
185 Int_t clusterizerType = -1;
186 Int_t eventType = -1;
188 clusterizerType = GetRecParam()->GetClusterizerFlag();
189 eventType = GetRecParam()->GetEventSpecie();
192 AliCDBEntry *entry = (AliCDBEntry*)
193 AliCDBManager::Instance()->Get("EMCAL/Calib/RecoParam");
194 //Get The reco param for the default event specie
196 AliEMCALRecParam *recParam = (AliEMCALRecParam*)((TObjArray *) entry->GetObject())->At(0);
197 if(recParam) clusterizerType = recParam->GetClusterizerFlag();
201 //Check if clusterizer previously set corresponds to what is needed for this event type
203 if(eventType!=AliRecoParam::kCalib){
204 //printf("ReCreate clusterizer? Clusterizer set <%d>, Clusterizer in use <%s>\n",
205 // clusterizerType, fgClusterizer->Version());
207 if (clusterizerType == AliEMCALRecParam::kClusterizerv1 && !strcmp(fgClusterizer->Version(),"clu-v1")) return;
209 else if(clusterizerType == AliEMCALRecParam::kClusterizerNxN && !strcmp(fgClusterizer->Version(),"clu-NxN")) return;
211 //Need to create new clusterizer, the one set previously is not the correct one
212 delete fgClusterizer;
217 if (clusterizerType == AliEMCALRecParam::kClusterizerv1)
219 fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData);
223 fgClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData,fPedestalData);
227 //____________________________________________________________________________
228 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
230 // method called by AliReconstruction;
231 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
232 // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by
233 // the global tracking.
234 // Works on the current event.
236 AliCodeTimerAuto("",0)
238 //Get input digits and put them in fgDigitsArr, clear the list before
239 ReadDigitsArrayFromTree(digitsTree);
243 fgClusterizer->InitParameters();
244 fgClusterizer->SetOutput(clustersTree);
246 //Skip clusterization of LED events
247 if (GetRecParam()->GetEventSpecie()!=AliRecoParam::kCalib){
249 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
251 fgClusterizer->SetInput(digitsTree);
253 //fgClusterizer->Digits2Clusters("deb all") ; //For debugging
254 fgClusterizer->Digits2Clusters("");
256 fgClusterizer->Clear();
258 }//digits array exists and has somethind
261 clustersTree->Fill();
263 // Deleting the recpoints at the end of the reconstruction call
264 fgClusterizer->DeleteRecPoints();
267 //____________________________________________________________________________
268 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
271 // Conversion from raw data to
273 // Works on a single-event basis
277 fTriggerData->SetMode(1);
279 if(fgDigitsArr) fgDigitsArr->Clear("C");
281 TClonesArray *digitsTrg = new TClonesArray("AliEMCALTriggerRawDigit", 32 * 96);
283 Int_t bufsize = 32000;
284 digitsTree->Branch("EMCAL", &fgDigitsArr, bufsize);
285 digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
287 //Skip calibration events do the rest
288 Bool_t doFit = kTRUE;
289 if ( !(GetRecParam()->FitLEDEvents()) && GetRecParam()->GetEventSpecie()==AliRecoParam::kCalib) doFit = kFALSE;
291 //must be done here because, in constructor, option is not yet known
292 fgRawUtils->SetOption(GetOption());
294 // fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
296 // fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
297 // fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
298 fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
299 fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
300 fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
301 fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
302 fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
304 //fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
305 //fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
307 // fgRawUtils->SetTimeMin(-99999 );
308 // fgRawUtils->SetTimeMax( 99999 );
310 fgRawUtils->Raw2Digits(rawReader,fgDigitsArr,fPedestalData,digitsTrg,fTriggerData);
312 }//skip calibration event
314 AliDebug(1," Calibration Event, skip!");
324 //____________________________________________________________________________
325 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
326 AliESDEvent* esd) const
328 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
330 // Works on the current event
331 // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
334 //########################################
336 //########################################
338 static int saveOnce = 0;
340 Int_t v0M[2] = {0, 0};
342 AliESDVZERO* esdV0 = esd->GetVZEROData();
346 for (Int_t i = 0; i < 32; i++)
348 v0M[0] += (Int_t)esdV0->GetAdcV0C(i);
349 v0M[1] += (Int_t)esdV0->GetAdcV0A(i);
354 AliWarning("Cannot retrieve V0 ESD! Run w/ null V0 charges");
357 if (fgTriggerDigits) fgTriggerDigits->Clear();
359 TBranch *branchtrg = digitsTree->GetBranch("EMTRG");
363 AliError("Can't get the branch with the EMCAL trigger digits!");
367 branchtrg->SetAddress(&fgTriggerDigits);
368 branchtrg->GetEntry(0);
370 // Note: fgTriggerProcessor reset done at the end of this method
371 fgTriggerProcessor->Digits2Trigger(fgTriggerDigits, v0M, fTriggerData);
374 AliESDCaloTrigger* trgESD = esd->GetCaloTrigger("EMCAL");
378 trgESD->Allocate(fgTriggerDigits->GetEntriesFast());
380 for (Int_t i = 0; i < fgTriggerDigits->GetEntriesFast(); i++)
382 AliEMCALTriggerRawDigit* rdig = (AliEMCALTriggerRawDigit*)fgTriggerDigits->At(i);
385 if (fGeom->GetPositionInEMCALFromAbsFastORIndex(rdig->GetId(), px, py))
387 Int_t a = -1, t = -1, times[10];
389 rdig->GetMaximum(a, t);
390 rdig->GetL0Times(times);
392 trgESD->Add(px, py, a, t, times, rdig->GetNL0Times(), rdig->GetL1TimeSum(), rdig->GetTriggerBits());
396 trgESD->SetL1Threshold(0, fTriggerData->GetL1GammaThreshold());
398 trgESD->SetL1Threshold(1, fTriggerData->GetL1JetThreshold() );
401 fTriggerData->GetL1V0(v0);
404 trgESD->SetL1FrameMask(fTriggerData->GetL1FrameMask());
406 if (!saveOnce && fTriggerData->GetL1DataDecoded())
409 fTriggerData->GetL1TriggerType(type);
411 esd->SetCaloTriggerType(type);
418 fTriggerData->Reset();
420 //########################################
421 //##############Fill CaloCells###############
422 //########################################
424 //Get input digits and put them in fgDigitsArr, clear the list before
425 ReadDigitsArrayFromTree(digitsTree);
427 Int_t nDigits = fgDigitsArr->GetEntries(), idignew = 0 ;
428 AliDebug(1,Form("%d digits",nDigits));
430 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
431 emcCells.CreateContainer(nDigits);
432 emcCells.SetType(AliVCaloCells::kEMCALCell);
434 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
435 const AliEMCALDigit * dig = (const AliEMCALDigit*)fgDigitsArr->At(idig);
436 if(dig->GetAmplitude() > 0 ){
437 energy = fgClusterizer->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
438 if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
439 emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());
444 emcCells.SetNumberOfCells(idignew);
447 //------------------------------------------------------------
448 //-----------------CLUSTERS-----------------------------
449 //------------------------------------------------------------
450 clustersTree->SetBranchStatus("*",0); //disable all branches
451 clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
452 if(fgClustersArr) fgClustersArr->Clear();
453 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
454 branch->SetAddress(&fgClustersArr);
456 //clustersTree->GetEvent(0);
458 Int_t nClusters = fgClustersArr->GetEntries(), nClustersNew=0;
459 AliDebug(1,Form("%d clusters",nClusters));
462 //########################################
463 //##############Fill CaloClusters#############
464 //########################################
465 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
466 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)fgClustersArr->At(iClust);
468 //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
469 // clust->Print(); //For debugging
470 // Get information from EMCAL reconstruction points
473 clust->GetGlobalPosition(gpos);
474 for (Int_t ixyz=0; ixyz<3; ixyz++)
475 xyz[ixyz] = gpos[ixyz];
477 clust->GetElipsAxis(elipAxis);
478 //Create digits lists
479 Int_t cellMult = clust->GetMultiplicity();
480 //TArrayS digiList(digitMult);
481 Float_t *amplFloat = clust->GetEnergiesList();
482 Int_t *digitInts = clust->GetAbsId();
483 TArrayS absIdList(cellMult);
484 TArrayD fracList(cellMult);
486 Int_t newCellMult = 0;
487 for (Int_t iCell=0; iCell<cellMult; iCell++) {
488 if (amplFloat[iCell] > 0) {
489 absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
491 if(emcCells.GetCellAmplitude(digitInts[iCell])>0 && GetRecParam()->GetUnfold()){
492 fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell]));//get cell calibration value
496 fracList[newCellMult] = 0;
502 absIdList.Set(newCellMult);
503 fracList.Set(newCellMult);
505 if(newCellMult > 0) { // accept cluster if it has some digit
508 Int_t parentMult = 0;
509 Int_t *parentList = clust->GetParents(parentMult);
510 // fills the ESDCaloCluster
511 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
512 ec->SetType(AliVCluster::kEMCALClusterv1);
513 ec->SetPosition(xyz);
514 ec->SetE(clust->GetEnergy());
516 //Distance to the nearest bad crystal
517 ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower());
519 ec->SetNCells(newCellMult);
520 //Change type of list from short to ushort
521 UShort_t *newAbsIdList = new UShort_t[newCellMult];
522 Double_t *newFracList = new Double_t[newCellMult];
523 for(Int_t i = 0; i < newCellMult ; i++) {
524 newAbsIdList[i]=absIdList[i];
525 newFracList[i] =fracList[i];
527 ec->SetCellsAbsId(newAbsIdList);
528 ec->SetCellsAmplitudeFraction(newFracList);
529 ec->SetDispersion(clust->GetDispersion());
530 ec->SetChi2(-1); //not yet implemented
531 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
532 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
533 ec->SetTOF(clust->GetTime()) ; //time-of-fligh
534 ec->SetNExMax(clust->GetNExMax()); //number of local maxima
537 TArrayI arrayParents(parentMult,parentList);
538 ec->AddLabels(arrayParents);
543 Int_t nTracks = esd->GetNumberOfTracks();
544 for (Int_t itrack = 0; itrack < nTracks; itrack++)
546 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
547 if(track->GetEMCALcluster()==iClust)
549 Double_t dEta=-999, dPhi=-999;
550 Bool_t isMatch = CalculateResidual(track, ec, dEta, dPhi);
553 // AliDebug(10, "Not good");
556 AliEMCALMatch *match = new AliEMCALMatch();
557 match->SetIndexT(itrack);
558 match->SetDistance(TMath::Sqrt(dEta*dEta+dPhi*dPhi));
559 match->SetdEta(dEta);
560 match->SetdPhi(dPhi);
561 fMatches->Add(match);
564 fMatches->Sort(kSortAscending); //Sort matched tracks from closest to furthest
565 Int_t nMatch = fMatches->GetEntries();
566 TArrayI arrayTrackMatched(nMatch);
567 for(Int_t imatch=0; imatch<nMatch; imatch++)
569 AliEMCALMatch *match = (AliEMCALMatch*)fMatches->At(imatch);
570 arrayTrackMatched[imatch] = match->GetIndexT();
573 ec->SetTrackDistance(match->GetdPhi(), match->GetdEta());
576 ec->AddTracksMatched(arrayTrackMatched);
578 //add the cluster to the esd object
579 esd->AddCaloCluster(ec);
582 delete [] newAbsIdList ;
583 delete [] newFracList ;
585 } // cycle on clusters
588 //Reset the index of matched cluster for tracks
589 //to the one in CaloCluster array
590 Int_t ncls = esd->GetNumberOfCaloClusters();
591 for(Int_t icl=0; icl<ncls; icl++)
593 AliESDCaloCluster *cluster = esd->GetCaloCluster(icl);
594 if(!cluster || !cluster->IsEMCAL()) continue;
595 TArrayI *trackIndex = cluster->GetTracksMatched();
596 for(Int_t itr=0; itr<trackIndex->GetSize(); itr++)
598 AliESDtrack *track = esd->GetTrack(trackIndex->At(itr));
599 track->SetEMCALcluster(cluster->GetID());
604 //Fill ESDCaloCluster with PID weights
605 AliEMCALPID *pid = new AliEMCALPID;
606 //pid->SetPrintInfo(kTRUE);
607 pid->SetReconstructor(kTRUE);
611 //Store EMCAL misalignment matrixes
612 FillMisalMatrixes(esd) ;
616 //==================================================================================
617 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
618 //Store EMCAL matrixes in ESD Header
620 //Check, if matrixes was already stored
621 for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
622 if(esd->GetEMCALMatrix(sm)!=0)
626 //Create and store matrixes
628 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
631 //Note, that owner of copied marixes will be header
632 const Int_t bufsize = 255;
634 TGeoHMatrix * m = 0x0;
635 for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
636 snprintf(path,bufsize,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
637 if(sm >= 10) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
639 if (gGeoManager->CheckPath(path)){
640 gGeoManager->cd(path);
641 m = gGeoManager->GetCurrentMatrix() ;
642 // printf("================================================= \n");
643 // printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
645 esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
646 // printf("================================================= \n");
649 esd->SetEMCALMatrix(NULL,sm) ;
654 //__________________________________________________________________________
655 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
657 // Read the digits from the input tree
658 // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
660 // Clear previous digits in the list
662 fgDigitsArr->Clear("C");
665 // It should not happen, but just in case ...
666 fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
669 // Read the digits from the input tree
670 TBranch *branch = digitsTree->GetBranch("EMCAL");
672 AliError("can't get the branch with the EMCAL digits !");
676 branch->SetAddress(&fgDigitsArr);
680 //==================================================================================
681 Bool_t AliEMCALReconstructor::CalculateResidual(AliESDtrack *track, AliESDCaloCluster *cluster, Double_t &dEta, Double_t &dPhi)const
684 // calculate the residual between track and cluster
687 // If the esdFriend is available, use the TPCOuter point as the starting point of extrapolation
688 // Otherwise use the TPCInner point
689 AliExternalTrackParam *trkParam = 0;
690 const AliESDfriendTrack* friendTrack = track->GetFriendTrack();
691 if(friendTrack && friendTrack->GetTPCOut())
692 trkParam = const_cast<AliExternalTrackParam*>(friendTrack->GetTPCOut());
694 trkParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
695 if(!trkParam) return kFALSE;
697 //Perform extrapolation
701 AliExternalTrackParam trkParamTmp (*trkParam);
702 cluster->GetPosition(clsPos);
703 TVector3 vec(clsPos[0],clsPos[1],clsPos[2]);
704 Double_t alpha = ((int)(vec.Phi()*TMath::RadToDeg()/20)+0.5)*20*TMath::DegToRad();
705 //Rotate to proper local coordinate
707 trkParamTmp.Rotate(alpha);
708 //extrapolation is done here
709 if(!AliTrackerBase::PropagateTrackToBxByBz(&trkParamTmp, vec.X(), track->GetMass(), GetRecParam()->GetExtrapolateStep(), kFALSE, 0.8, -1))
712 //Calculate the residuals
713 trkParamTmp.GetXYZ(trkPos);
715 TVector3 clsPosVec(clsPos[0],clsPos[1],clsPos[2]);
716 TVector3 trkPosVec(trkPos[0],trkPos[1],trkPos[2]);
718 Double_t clsPhi = clsPosVec.Phi();
719 if(clsPhi<0) clsPhi+=2*TMath::Pi();
720 Double_t trkPhi = trkPosVec.Phi();
721 if(trkPhi<0) trkPhi+=2*TMath::Pi();
723 dPhi = clsPhi-trkPhi;
724 dEta = clsPosVec.Eta()-trkPosVec.Eta();
730 //==================================================================================
732 AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch()
739 //default constructor
744 //==================================================================================
746 AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch(const AliEMCALMatch& copy)
748 fIndexT(copy.fIndexT),
749 fDistance(copy.fDistance),
757 //==================================================================================
759 Int_t AliEMCALReconstructor::AliEMCALMatch::Compare(const TObject *obj) const
762 // Compare wrt the residual
765 AliEMCALReconstructor::AliEMCALMatch *that = (AliEMCALReconstructor::AliEMCALMatch*)obj;
767 Double_t thisDist = fDistance;//fDistance;
768 Double_t thatDist = that->fDistance;//that->GetDistance();
770 if (thisDist > thatDist) return 1;
771 else if (thisDist < thatDist) return -1;