00001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00002
00003 #include "RecoParticleFlow/PFProducer/interface/PFAlgo.h"
00004 #include "RecoParticleFlow/PFProducer/interface/PFMuonAlgo.h"
00005 #include "RecoParticleFlow/PFProducer/interface/PFElectronAlgo.h"
00006 #include "RecoParticleFlow/PFProducer/interface/PFConversionAlgo.h"
00007 #include "RecoParticleFlow/PFProducer/interface/PFPhotonAlgo.h"
00008 #include "RecoParticleFlow/PFProducer/interface/PFElectronExtraEqual.h"
00009
00010 #include "RecoParticleFlow/PFClusterTools/interface/PFEnergyCalibration.h"
00011 #include "RecoParticleFlow/PFClusterTools/interface/PFEnergyCalibrationHF.h"
00012 #include "RecoParticleFlow/PFClusterTools/interface/PFSCEnergyCalibration.h"
00013
00014 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00015 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00016 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
00017 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
00018 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
00019 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00020 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
00021
00022 #include "DataFormats/VertexReco/interface/Vertex.h"
00023 #include "DataFormats/TrackReco/interface/Track.h"
00024 #include "DataFormats/MuonReco/interface/Muon.h"
00025
00026
00027 #include "DataFormats/Common/interface/OrphanHandle.h"
00028
00029 #include "DataFormats/Provenance/interface/ProductID.h"
00030 #include "DataFormats/Math/interface/LorentzVector.h"
00031
00032
00033 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00034 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00035
00036 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertex.h"
00037 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertexFwd.h"
00038
00039
00040 #include "Math/PxPyPzM4D.h"
00041 #include "Math/LorentzVector.h"
00042 #include "Math/DisplacementVector3D.h"
00043 #include "Math/SMatrix.h"
00044 #include "TDecompChol.h"
00045
00046 #include "boost/graph/adjacency_matrix.hpp"
00047 #include "boost/graph/graph_utility.hpp"
00048
00049
00050 using namespace std;
00051 using namespace reco;
00052 using namespace boost;
00053
00054 typedef std::list< reco::PFBlockRef >::iterator IBR;
00055
00056
00057
00058 PFAlgo::PFAlgo()
00059 : pfCandidates_( new PFCandidateCollection),
00060 nSigmaECAL_(0),
00061 nSigmaHCAL_(1),
00062 algo_(1),
00063 debug_(false),
00064 pfele_(0),
00065 pfpho_(0),
00066 useVertices_(false)
00067 {}
00068
00069 PFAlgo::~PFAlgo() {
00070 if (usePFElectrons_) delete pfele_;
00071 if (usePFConversions_) delete pfConversion_;
00072 if (usePFPhotons_) delete pfpho_;
00073 }
00074
00075
00076 void
00077 PFAlgo::setParameters(double nSigmaECAL,
00078 double nSigmaHCAL,
00079 const boost::shared_ptr<PFEnergyCalibration>& calibration,
00080 const boost::shared_ptr<PFEnergyCalibrationHF>& thepfEnergyCalibrationHF) {
00081
00082 nSigmaECAL_ = nSigmaECAL;
00083 nSigmaHCAL_ = nSigmaHCAL;
00084
00085 calibration_ = calibration;
00086 thepfEnergyCalibrationHF_ = thepfEnergyCalibrationHF;
00087
00088 }
00089
00090
00091 void
00092 PFAlgo::setPFEleParameters(double mvaEleCut,
00093 string mvaWeightFileEleID,
00094 bool usePFElectrons,
00095 const boost::shared_ptr<PFSCEnergyCalibration>& thePFSCEnergyCalibration,
00096 const boost::shared_ptr<PFEnergyCalibration>& thePFEnergyCalibration,
00097 double sumEtEcalIsoForEgammaSC_barrel,
00098 double sumEtEcalIsoForEgammaSC_endcap,
00099 double coneEcalIsoForEgammaSC,
00100 double sumPtTrackIsoForEgammaSC_barrel,
00101 double sumPtTrackIsoForEgammaSC_endcap,
00102 unsigned int nTrackIsoForEgammaSC,
00103 double coneTrackIsoForEgammaSC,
00104 bool applyCrackCorrections,
00105 bool usePFSCEleCalib,
00106 bool useEGElectrons,
00107 bool useEGammaSupercluster) {
00108
00109 mvaEleCut_ = mvaEleCut;
00110 usePFElectrons_ = usePFElectrons;
00111 applyCrackCorrectionsElectrons_ = applyCrackCorrections;
00112 usePFSCEleCalib_ = usePFSCEleCalib;
00113 thePFSCEnergyCalibration_ = thePFSCEnergyCalibration;
00114 useEGElectrons_ = useEGElectrons;
00115 useEGammaSupercluster_ = useEGammaSupercluster;
00116 sumEtEcalIsoForEgammaSC_barrel_ = sumEtEcalIsoForEgammaSC_barrel;
00117 sumEtEcalIsoForEgammaSC_endcap_ = sumEtEcalIsoForEgammaSC_endcap;
00118 coneEcalIsoForEgammaSC_ = coneEcalIsoForEgammaSC;
00119 sumPtTrackIsoForEgammaSC_barrel_ = sumPtTrackIsoForEgammaSC_barrel;
00120 sumPtTrackIsoForEgammaSC_endcap_ = sumPtTrackIsoForEgammaSC_endcap;
00121 coneTrackIsoForEgammaSC_ = coneTrackIsoForEgammaSC;
00122 nTrackIsoForEgammaSC_ = nTrackIsoForEgammaSC;
00123
00124
00125 if(!usePFElectrons_) return;
00126 mvaWeightFileEleID_ = mvaWeightFileEleID;
00127 FILE * fileEleID = fopen(mvaWeightFileEleID_.c_str(), "r");
00128 if (fileEleID) {
00129 fclose(fileEleID);
00130 }
00131 else {
00132 string err = "PFAlgo: cannot open weight file '";
00133 err += mvaWeightFileEleID;
00134 err += "'";
00135 throw invalid_argument( err );
00136 }
00137 pfele_= new PFElectronAlgo(mvaEleCut_,mvaWeightFileEleID_,
00138 thePFSCEnergyCalibration_,
00139 thePFEnergyCalibration,
00140 applyCrackCorrectionsElectrons_,
00141 usePFSCEleCalib_,
00142 useEGElectrons_,
00143 useEGammaSupercluster_,
00144 sumEtEcalIsoForEgammaSC_barrel_,
00145 sumEtEcalIsoForEgammaSC_endcap_,
00146 coneEcalIsoForEgammaSC_,
00147 sumPtTrackIsoForEgammaSC_barrel_,
00148 sumPtTrackIsoForEgammaSC_endcap_,
00149 nTrackIsoForEgammaSC_,
00150 coneTrackIsoForEgammaSC_);
00151 }
00152
00153 void
00154 PFAlgo::setPFPhotonParameters(bool usePFPhotons,
00155 std::string mvaWeightFileConvID,
00156 double mvaConvCut,
00157 const boost::shared_ptr<PFEnergyCalibration>& thePFEnergyCalibration) {
00158
00159 usePFPhotons_ = usePFPhotons;
00160
00161
00162 reco::Vertex dummy;
00163 const reco::Vertex* pv=&dummy;
00164 if(useVertices_)
00165 {
00166 dummy = primaryVertex_;
00167 }
00168 else {
00169 reco::Vertex::Error e;
00170 e(0, 0) = 0.0015 * 0.0015;
00171 e(1, 1) = 0.0015 * 0.0015;
00172 e(2, 2) = 15. * 15.;
00173 reco::Vertex::Point p(0, 0, 0);
00174 dummy = reco::Vertex(p, e, 0, 0, 0);
00175 }
00176
00177 if(! usePFPhotons_) return;
00178 FILE * filePhotonConvID = fopen(mvaWeightFileConvID.c_str(), "r");
00179 if (filePhotonConvID) {
00180 fclose(filePhotonConvID);
00181 }
00182 else {
00183 string err = "PFAlgo: cannot open weight file '";
00184 err += mvaWeightFileConvID;
00185 err += "'";
00186 throw invalid_argument( err );
00187 }
00188
00189
00190 pfpho_ = new PFPhotonAlgo(mvaWeightFileConvID,
00191 mvaConvCut,
00192 *pv,
00193 thePFEnergyCalibration);
00194 return;
00195 }
00196
00197 void
00198 PFAlgo::setPFMuonAndFakeParameters(std::vector<double> muonHCAL,
00199 std::vector<double> muonECAL,
00200 double nSigmaTRACK,
00201 double ptError,
00202 std::vector<double> factors45,
00203 bool usePFMuonMomAssign)
00204 {
00205 muonHCAL_ = muonHCAL;
00206 muonECAL_ = muonECAL;
00207 nSigmaTRACK_ = nSigmaTRACK;
00208 ptError_ = ptError;
00209 factors45_ = factors45;
00210 usePFMuonMomAssign_ = usePFMuonMomAssign;
00211 }
00212
00213 void
00214 PFAlgo::setPostHFCleaningParameters(bool postHFCleaning,
00215 double minHFCleaningPt,
00216 double minSignificance,
00217 double maxSignificance,
00218 double minSignificanceReduction,
00219 double maxDeltaPhiPt,
00220 double minDeltaMet) {
00221 postHFCleaning_ = postHFCleaning;
00222 minHFCleaningPt_ = minHFCleaningPt;
00223 minSignificance_ = minSignificance;
00224 maxSignificance_ = maxSignificance;
00225 minSignificanceReduction_= minSignificanceReduction;
00226 maxDeltaPhiPt_ = maxDeltaPhiPt;
00227 minDeltaMet_ = minDeltaMet;
00228 }
00229
00230 void
00231 PFAlgo::setDisplacedVerticesParameters(bool rejectTracks_Bad,
00232 bool rejectTracks_Step45,
00233 bool usePFNuclearInteractions,
00234 bool usePFConversions,
00235 bool usePFDecays,
00236 double dptRel_DispVtx){
00237
00238 rejectTracks_Bad_ = rejectTracks_Bad;
00239 rejectTracks_Step45_ = rejectTracks_Step45;
00240 usePFNuclearInteractions_ = usePFNuclearInteractions;
00241 usePFConversions_ = usePFConversions;
00242 pfConversion_ = new PFConversionAlgo();
00243 usePFDecays_ = usePFDecays;
00244 dptRel_DispVtx_ = dptRel_DispVtx;
00245
00246 }
00247
00248
00249 void
00250 PFAlgo::setPFVertexParameters(bool useVertex,
00251 const reco::VertexCollection& primaryVertices) {
00252 useVertices_ = useVertex;
00253
00254 bool primaryVertexFound = false;
00255 for (unsigned short i=0 ;i<primaryVertices.size();++i)
00256 {
00257 if(primaryVertices[i].isValid()&&(!primaryVertices[i].isFake()))
00258 {
00259 primaryVertex_ = primaryVertices[i];
00260 primaryVertexFound = true;
00261 break;
00262 }
00263 }
00264
00265 useVertices_ = useVertex && primaryVertexFound;
00266
00267 }
00268
00269
00270 void PFAlgo::reconstructParticles( const reco::PFBlockHandle& blockHandle ) {
00271
00272 blockHandle_ = blockHandle;
00273 reconstructParticles( *blockHandle_ );
00274 }
00275
00276
00277
00278 void PFAlgo::reconstructParticles( const reco::PFBlockCollection& blocks ) {
00279
00280
00281 if(pfCandidates_.get() )
00282 pfCandidates_->clear();
00283 else
00284 pfCandidates_.reset( new reco::PFCandidateCollection );
00285
00286 if(pfElectronCandidates_.get() )
00287 pfElectronCandidates_->clear();
00288 else
00289 pfElectronCandidates_.reset( new reco::PFCandidateCollection);
00290
00291
00292 if( pfPhotonCandidates_.get() )
00293 pfPhotonCandidates_->clear();
00294 else
00295 pfPhotonCandidates_.reset( new reco::PFCandidateCollection);
00296
00297 if(pfCleanedCandidates_.get() )
00298 pfCleanedCandidates_->clear();
00299 else
00300 pfCleanedCandidates_.reset( new reco::PFCandidateCollection );
00301
00302 if(pfCosmicsMuonCleanedCandidates_.get() )
00303 pfCosmicsMuonCleanedCandidates_->clear();
00304 else
00305 pfCosmicsMuonCleanedCandidates_.reset( new reco::PFCandidateCollection );
00306
00307 if(pfCleanedTrackerAndGlobalMuonCandidates_.get() )
00308 pfCleanedTrackerAndGlobalMuonCandidates_->clear();
00309 else
00310 pfCleanedTrackerAndGlobalMuonCandidates_.reset( new reco::PFCandidateCollection );
00311
00312 if(pfFakeMuonCleanedCandidates_.get() )
00313 pfFakeMuonCleanedCandidates_->clear();
00314 else
00315 pfFakeMuonCleanedCandidates_.reset( new reco::PFCandidateCollection );
00316
00317 if(pfPunchThroughMuonCleanedCandidates_.get() )
00318 pfPunchThroughMuonCleanedCandidates_->clear();
00319 else
00320 pfPunchThroughMuonCleanedCandidates_.reset( new reco::PFCandidateCollection );
00321
00322 if(pfPunchThroughHadronCleanedCandidates_.get() )
00323 pfPunchThroughHadronCleanedCandidates_->clear();
00324 else
00325 pfPunchThroughHadronCleanedCandidates_.reset( new reco::PFCandidateCollection );
00326
00327 if(pfAddedMuonCandidates_.get() )
00328 pfAddedMuonCandidates_->clear();
00329 else
00330 pfAddedMuonCandidates_.reset( new reco::PFCandidateCollection );
00331
00332
00333 pfElectronExtra_.clear();
00334
00335 if( debug_ ) {
00336 cout<<"*********************************************************"<<endl;
00337 cout<<"***** Particle flow algorithm *****"<<endl;
00338 cout<<"*********************************************************"<<endl;
00339 }
00340
00341
00342 std::list< reco::PFBlockRef > hcalBlockRefs;
00343 std::list< reco::PFBlockRef > ecalBlockRefs;
00344 std::list< reco::PFBlockRef > otherBlockRefs;
00345
00346 for( unsigned i=0; i<blocks.size(); ++i ) {
00347
00348 reco::PFBlockRef blockref = createBlockRef( blocks, i);
00349
00350 const reco::PFBlock& block = *blockref;
00351 const edm::OwnVector< reco::PFBlockElement >&
00352 elements = block.elements();
00353
00354 bool singleEcalOrHcal = false;
00355 if( elements.size() == 1 ){
00356 if( elements[0].type() == reco::PFBlockElement::ECAL ){
00357 ecalBlockRefs.push_back( blockref );
00358 singleEcalOrHcal = true;
00359 }
00360 if( elements[0].type() == reco::PFBlockElement::HCAL ){
00361 hcalBlockRefs.push_back( blockref );
00362 singleEcalOrHcal = true;
00363 }
00364 }
00365
00366 if(!singleEcalOrHcal) {
00367 otherBlockRefs.push_back( blockref );
00368 }
00369 }
00370
00371 if( debug_ ){
00372 cout<<"# Ecal blocks: "<<ecalBlockRefs.size()
00373 <<", # Hcal blocks: "<<hcalBlockRefs.size()
00374 <<", # Other blocks: "<<otherBlockRefs.size()<<endl;
00375 }
00376
00377
00378
00379
00380
00381 for( IBR io = otherBlockRefs.begin(); io!=otherBlockRefs.end(); ++io) {
00382 processBlock( *io, hcalBlockRefs, ecalBlockRefs );
00383 }
00384
00385 std::list< reco::PFBlockRef > empty;
00386
00387
00388 for( IBR ih = hcalBlockRefs.begin(); ih!=hcalBlockRefs.end(); ++ih) {
00389 processBlock( *ih, empty, empty );
00390 }
00391
00392
00393 for( IBR ie = ecalBlockRefs.begin(); ie!=ecalBlockRefs.end(); ++ie) {
00394 processBlock( *ie, empty, empty );
00395 }
00396
00397
00398 postCleaning();
00399 }
00400
00401
00402 void PFAlgo::processBlock( const reco::PFBlockRef& blockref,
00403 std::list<reco::PFBlockRef>& hcalBlockRefs,
00404 std::list<reco::PFBlockRef>& ecalBlockRefs ) {
00405
00406
00407 assert(!blockref.isNull() );
00408 const reco::PFBlock& block = *blockref;
00409
00410 typedef std::multimap<double, unsigned>::iterator IE;
00411 typedef std::multimap<double, std::pair<unsigned,math::XYZVector> >::iterator IS;
00412 typedef std::multimap<double, std::pair<unsigned,bool> >::iterator IT;
00413 typedef std::multimap< unsigned, std::pair<double, unsigned> >::iterator II;
00414
00415 if(debug_) {
00416 cout<<"#########################################################"<<endl;
00417 cout<<"##### Process Block: #####"<<endl;
00418 cout<<"#########################################################"<<endl;
00419 cout<<block<<endl;
00420 }
00421
00422 const edm::OwnVector< reco::PFBlockElement >& elements = block.elements();
00423
00424
00425 PFBlock::LinkData linkData = block.linkData();
00426
00427
00428 vector<bool> active( elements.size(), true );
00429
00430
00431
00432
00433 if (usePFElectrons_) {
00434 if (pfele_->isElectronValidCandidate(blockref,active)){
00435
00436 const std::vector<reco::PFCandidate> & PFElectCandidates_(pfele_->getElectronCandidates());
00437
00438
00439
00440
00441 for ( std::vector<reco::PFCandidate>::const_iterator ec = PFElectCandidates_.begin();
00442 ec != PFElectCandidates_.end(); ++ec ) {
00443 pfCandidates_->push_back(*ec);
00444
00445 }
00446
00447
00448 }
00449 pfElectronCandidates_->insert(pfElectronCandidates_->end(),
00450 pfele_->getAllElectronCandidates().begin(),
00451 pfele_->getAllElectronCandidates().end());
00452
00453
00454 pfElectronExtra_.insert(pfElectronExtra_.end(),
00455 pfele_->getElectronExtra().begin(),
00456 pfele_->getElectronExtra().end());
00457 }
00458
00459
00460
00461 if (usePFConversions_ && false) {
00462 if (pfConversion_->isConversionValidCandidate(blockref, active )){
00463
00464 std::vector<reco::PFCandidate> PFConversionCandidates = pfConversion_->conversionCandidates();
00465
00466 for ( std::vector<reco::PFCandidate>::iterator iConv = PFConversionCandidates.begin();
00467 iConv != PFConversionCandidates.end(); ++iConv ) {
00468 pfCandidates_->push_back(*iConv);
00469 }
00470 for(unsigned iEle=0; iEle<elements.size(); iEle++) {
00471
00472 }
00473 }
00474 }
00475
00476 if( usePFPhotons_ ) {
00477
00478 if(debug_)
00479 cout<<endl<<"--------------- entering PFPhotonAlgo ----------------"<<endl;
00480
00481 if ( pfpho_->isPhotonValidCandidate(blockref,
00482 active,
00483 pfPhotonCandidates_,
00484 pfElectronCandidates_
00485 ) ) {
00486 if(debug_)
00487 std::cout<< " In this PFBlock we found "<<pfPhotonCandidates_->size()<<" Photon Candidates."<<std::endl;
00488
00489
00490
00491
00492
00493 PFCandidateCollection::const_iterator cand = pfPhotonCandidates_->begin();
00494 for( ; cand != pfPhotonCandidates_->end(); ++cand)
00495 pfCandidates_->push_back(*cand);
00496
00497 }
00498 pfPhotonCandidates_->clear();
00499 }
00500
00501 if(debug_)
00502 cout<<endl<<"--------------- loop 1 ------------------"<<endl;
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517
00518
00519
00520
00521
00522
00523
00524
00525
00526
00527
00528
00529 vector<unsigned> hcalIs;
00530 vector<unsigned> ecalIs;
00531 vector<unsigned> trackIs;
00532 vector<unsigned> ps1Is;
00533 vector<unsigned> ps2Is;
00534
00535 vector<unsigned> hfEmIs;
00536 vector<unsigned> hfHadIs;
00537
00538 for(unsigned iEle=0; iEle<elements.size(); iEle++) {
00539
00540 PFBlockElement::Type type = elements[iEle].type();
00541
00542 if(debug_ && type != PFBlockElement::BREM ) cout<<endl<<elements[iEle];
00543
00544 switch( type ) {
00545 case PFBlockElement::TRACK:
00546 if ( active[iEle] ) {
00547 trackIs.push_back( iEle );
00548 if(debug_) cout<<"TRACK, stored index, continue"<<endl;
00549 }
00550 break;
00551 case PFBlockElement::ECAL:
00552 if ( active[iEle] ) {
00553 ecalIs.push_back( iEle );
00554 if(debug_) cout<<"ECAL, stored index, continue"<<endl;
00555 }
00556 continue;
00557 case PFBlockElement::HCAL:
00558 if ( active[iEle] ) {
00559 hcalIs.push_back( iEle );
00560 if(debug_) cout<<"HCAL, stored index, continue"<<endl;
00561 }
00562 continue;
00563 case PFBlockElement::HFEM:
00564 if ( active[iEle] ) {
00565 hfEmIs.push_back( iEle );
00566 if(debug_) cout<<"HFEM, stored index, continue"<<endl;
00567 }
00568 continue;
00569 case PFBlockElement::HFHAD:
00570 if ( active[iEle] ) {
00571 hfHadIs.push_back( iEle );
00572 if(debug_) cout<<"HFHAD, stored index, continue"<<endl;
00573 }
00574 continue;
00575
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585
00586
00587
00588
00589
00590 default:
00591 continue;
00592 }
00593
00594
00595 unsigned iTrack = iEle;
00596 reco::MuonRef muonRef = elements[iTrack].muonRef();
00597
00598
00599 if (active[iTrack] && isFromSecInt(elements[iEle], "primary")){
00600 bool isPrimaryTrack = elements[iEle].displacedVertexRef(PFBlockElement::T_TO_DISP)->displacedVertexRef()->isTherePrimaryTracks();
00601 if (isPrimaryTrack) {
00602 if (debug_) cout << "Primary Track reconstructed alone" << endl;
00603 unsigned tmpi = reconstructTrack(elements[iEle]);
00604 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEle );
00605 active[iTrack] = false;
00606 }
00607 }
00608
00609
00610 if(debug_) {
00611 if ( !active[iTrack] )
00612 cout << "Already used by electrons, muons, conversions" << endl;
00613 }
00614
00615
00616
00617 if ( ! active[iTrack] ) continue;
00618
00619 reco::TrackRef trackRef = elements[iTrack].trackRef();
00620 assert( !trackRef.isNull() );
00621
00622
00623
00624
00625
00626
00627
00628
00629 std::multimap<double, unsigned> ecalElems;
00630 block.associatedElements( iTrack, linkData,
00631 ecalElems ,
00632 reco::PFBlockElement::ECAL,
00633 reco::PFBlock::LINKTEST_ALL );
00634
00635 std::multimap<double, unsigned> hcalElems;
00636 block.associatedElements( iTrack, linkData,
00637 hcalElems,
00638 reco::PFBlockElement::HCAL,
00639 reco::PFBlock::LINKTEST_ALL );
00640
00641
00642
00643
00644 if ( hcalElems.empty() && !ecalElems.empty() ) {
00645
00646 unsigned ntt = 1;
00647 unsigned index = ecalElems.begin()->second;
00648 std::multimap<double, unsigned> sortedTracks;
00649 block.associatedElements( index, linkData,
00650 sortedTracks,
00651 reco::PFBlockElement::TRACK,
00652 reco::PFBlock::LINKTEST_ALL );
00653
00654
00655
00656 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
00657 unsigned jTrack = ie->second;
00658
00659
00660 if ( !active[jTrack] ) continue;
00661
00662
00663
00664 if ( jTrack == iTrack ) continue;
00665
00666
00667
00668
00669 std::multimap<double, unsigned> sortedECAL;
00670 block.associatedElements( jTrack, linkData,
00671 sortedECAL,
00672 reco::PFBlockElement::ECAL,
00673 reco::PFBlock::LINKTEST_ALL );
00674 if ( sortedECAL.begin()->second != index ) continue;
00675
00676
00677
00678
00679 std::multimap<double, unsigned> sortedHCAL;
00680 block.associatedElements( jTrack, linkData,
00681 sortedHCAL,
00682 reco::PFBlockElement::HCAL,
00683 reco::PFBlock::LINKTEST_ALL );
00684 if ( !sortedHCAL.size() ) continue;
00685
00686 ntt++;
00687
00688
00689 block.setLink( iTrack,
00690 sortedHCAL.begin()->second,
00691 sortedECAL.begin()->first,
00692 linkData,
00693 PFBlock::LINKTEST_RECHIT );
00694
00695 }
00696
00697
00698 block.associatedElements( iTrack, linkData,
00699 hcalElems,
00700 reco::PFBlockElement::HCAL,
00701 reco::PFBlock::LINKTEST_ALL );
00702
00703 if ( debug_ && hcalElems.size() )
00704 std::cout << "Track linked back to HCAL due to ECAL sharing with other tracks" << std::endl;
00705 }
00706
00707
00708
00709
00710
00711
00712
00713
00714 std::multimap<double,unsigned> gsfElems;
00715 if (usePFElectrons_ == false) {
00716 block.associatedElements( iTrack, linkData,
00717 gsfElems ,
00718 reco::PFBlockElement::GSF );
00719 }
00720
00721
00722 if(hcalElems.empty() && debug_) {
00723 cout<<"no hcal element connected to track "<<iTrack<<endl;
00724 }
00725
00726
00727
00728
00729 bool hcalFound = false;
00730
00731 if(debug_)
00732 cout<<"now looping on elements associated to the track"<<endl;
00733
00734
00735
00736 for(IE ie = ecalElems.begin(); ie != ecalElems.end(); ++ie ) {
00737
00738 unsigned index = ie->second;
00739
00740 PFBlockElement::Type type = elements[index].type();
00741 if(debug_) {
00742 double dist = ie->first;
00743 cout<<"\telement "<<elements[index]<<" linked with distance = "<< dist <<endl;
00744 if ( ! active[index] ) cout << "This ECAL is already used - skip it" << endl;
00745 }
00746 assert( type == PFBlockElement::ECAL );
00747
00748
00749 if ( ! active[index] ) continue;
00750
00751
00752
00753
00754
00755 if( !hcalElems.empty() && debug_)
00756 cout<<"\t\tat least one hcal element connected to the track."
00757 <<" Sparing Ecal cluster for the hcal loop"<<endl;
00758
00759 }
00760
00761
00762
00763
00764
00765 if( hcalElems.empty() ) {
00766
00767 if ( debug_ )
00768 std::cout << "Now deals with tracks linked to no HCAL clusters" << std::endl;
00769
00770
00771
00772
00773 double trackMomentum = elements[iTrack].trackRef()->p();
00774 if ( debug_ )
00775 std::cout << elements[iTrack] << std::endl;
00776
00777
00778 bool thisIsAMuon = PFMuonAlgo::isMuon(elements[iTrack]);
00779 if ( thisIsAMuon ) trackMomentum = 0.;
00780
00781
00782 bool rejectFake = false;
00783 if ( !thisIsAMuon && elements[iTrack].trackRef()->ptError() > ptError_ ) {
00784
00785 double deficit = trackMomentum;
00786 double resol = neutralHadronEnergyResolution(trackMomentum,
00787 elements[iTrack].trackRef()->eta());
00788 resol *= trackMomentum;
00789
00790 if ( !ecalElems.empty() ) {
00791 unsigned thisEcal = ecalElems.begin()->second;
00792 reco::PFClusterRef clusterRef = elements[thisEcal].clusterRef();
00793 deficit -= clusterRef->energy();
00794 resol = neutralHadronEnergyResolution(trackMomentum,
00795 clusterRef->positionREP().Eta());
00796 resol *= trackMomentum;
00797 }
00798
00799 bool isPrimary = isFromSecInt(elements[iTrack], "primary");
00800
00801 if ( deficit > nSigmaTRACK_*resol && !isPrimary) {
00802 rejectFake = true;
00803 active[iTrack] = false;
00804 if ( debug_ )
00805 std::cout << elements[iTrack] << std::endl
00806 << "is probably a fake (1) --> lock the track"
00807 << std::endl;
00808 }
00809 }
00810
00811 if ( rejectFake ) continue;
00812
00813
00814
00815
00816 std::vector<unsigned> tmpi;
00817 std::vector<unsigned> kTrack;
00818
00819
00820 double Dpt = trackRef->ptError();
00821
00822 if ( rejectTracks_Step45_ && ecalElems.empty() &&
00823 trackMomentum > 30. && Dpt > 0.5 &&
00824 ( trackRef->algo() == TrackBase::iter3 ||
00825 trackRef->algo() == TrackBase::iter4 ||
00826 trackRef->algo() == TrackBase::iter5 ) ) {
00827
00828
00829 double dptRel = Dpt/trackRef->pt()*100;
00830 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
00831
00832 if ( isPrimaryOrSecondary && dptRel < dptRel_DispVtx_) continue;
00833 unsigned nHits = elements[iTrack].trackRef()->hitPattern().trackerLayersWithMeasurement();
00834 unsigned int NLostHit = trackRef->hitPattern().trackerLayersWithoutMeasurement();
00835
00836 if ( debug_ )
00837 std::cout << "A track (algo = " << trackRef->algo() << ") with momentum " << trackMomentum
00838 << " / " << elements[iTrack].trackRef()->pt() << " +/- " << Dpt
00839 << " / " << elements[iTrack].trackRef()->eta()
00840 << " without any link to ECAL/HCAL and with " << nHits << " (" << NLostHit
00841 << ") hits (lost hits) has been cleaned" << std::endl;
00842 active[iTrack] = false;
00843 continue;
00844 }
00845
00846
00847 tmpi.push_back(reconstructTrack( elements[iTrack]));
00848 kTrack.push_back(iTrack);
00849 active[iTrack] = false;
00850
00851
00852 if ( ecalElems.empty() ) {
00853 (*pfCandidates_)[tmpi[0]].setEcalEnergy( 0., 0. );
00854 (*pfCandidates_)[tmpi[0]].setHcalEnergy( 0., 0. );
00855 (*pfCandidates_)[tmpi[0]].setPs1Energy( 0 );
00856 (*pfCandidates_)[tmpi[0]].setPs2Energy( 0 );
00857 (*pfCandidates_)[tmpi[0]].addElementInBlock( blockref, kTrack[0] );
00858 continue;
00859 }
00860
00861
00862 unsigned thisEcal = ecalElems.begin()->second;
00863 reco::PFClusterRef clusterRef = elements[thisEcal].clusterRef();
00864 if ( debug_ ) std::cout << " is associated to " << elements[thisEcal] << std::endl;
00865
00866
00867 if ( thisIsAMuon ) {
00868 (*pfCandidates_)[tmpi[0]].setEcalEnergy( clusterRef->energy(),
00869 std::min(clusterRef->energy(), muonECAL_[0]) );
00870 (*pfCandidates_)[tmpi[0]].setHcalEnergy( 0., 0. );
00871 (*pfCandidates_)[tmpi[0]].setPs1Energy( 0 );
00872 (*pfCandidates_)[tmpi[0]].setPs2Energy( 0 );
00873 (*pfCandidates_)[tmpi[0]].addElementInBlock( blockref, kTrack[0] );
00874 }
00875
00876 double slopeEcal = 1.;
00877 bool connectedToEcal = false;
00878 unsigned iEcal = 99999;
00879 double calibEcal = 0.;
00880 double calibHcal = 0.;
00881 double totalEcal = thisIsAMuon ? -muonECAL_[0] : 0.;
00882 double totalHcal = 0.;
00883
00884
00885 std::multimap<double, unsigned> sortedTracks;
00886 block.associatedElements( thisEcal, linkData,
00887 sortedTracks,
00888 reco::PFBlockElement::TRACK,
00889 reco::PFBlock::LINKTEST_ALL );
00890
00891 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
00892 unsigned jTrack = ie->second;
00893
00894
00895 if ( !active[jTrack] ) continue;
00896
00897
00898 if ( jTrack == iTrack ) continue;
00899
00900
00901
00902
00903 std::multimap<double, unsigned> sortedECAL;
00904 block.associatedElements( jTrack, linkData,
00905 sortedECAL,
00906 reco::PFBlockElement::ECAL,
00907 reco::PFBlock::LINKTEST_ALL );
00908 if ( sortedECAL.begin()->second != thisEcal ) continue;
00909
00910
00911 bool thatIsAMuon = PFMuonAlgo::isMuon(elements[jTrack]);
00912
00913
00914 bool rejectFake = false;
00915 reco::TrackRef trackRef = elements[jTrack].trackRef();
00916 if ( !thatIsAMuon && trackRef->ptError() > ptError_) {
00917 double deficit = trackMomentum + trackRef->p() - clusterRef->energy();
00918 double resol = nSigmaTRACK_*neutralHadronEnergyResolution(trackMomentum+trackRef->p(),
00919 clusterRef->positionREP().Eta());
00920 resol *= (trackMomentum+trackRef->p());
00921 if ( deficit > nSigmaTRACK_*resol ) {
00922 rejectFake = true;
00923 kTrack.push_back(jTrack);
00924 active[jTrack] = false;
00925 if ( debug_ )
00926 std::cout << elements[jTrack] << std::endl
00927 << "is probably a fake (2) --> lock the track"
00928 << std::endl;
00929 }
00930 }
00931 if ( rejectFake ) continue;
00932
00933
00934
00935
00936
00937 if ( !thatIsAMuon ) {
00938 if ( debug_ )
00939 std::cout << "Track momentum increased from " << trackMomentum << " GeV ";
00940 trackMomentum += trackRef->p();
00941 if ( debug_ ) {
00942 std::cout << "to " << trackMomentum << " GeV." << std::endl;
00943 std::cout << "with " << elements[jTrack] << std::endl;
00944 }
00945 } else {
00946 totalEcal -= muonECAL_[0];
00947 totalEcal = std::max(totalEcal, 0.);
00948 }
00949
00950
00951 tmpi.push_back(reconstructTrack( elements[jTrack] ));
00952 kTrack.push_back(jTrack);
00953 active[jTrack] = false;
00954
00955 if ( thatIsAMuon ) {
00956 (*pfCandidates_)[tmpi.back()].setEcalEnergy(clusterRef->energy(),
00957 std::min(clusterRef->energy(),muonECAL_[0]));
00958 (*pfCandidates_)[tmpi.back()].setHcalEnergy( 0., 0. );
00959 (*pfCandidates_)[tmpi.back()].setPs1Energy( 0 );
00960 (*pfCandidates_)[tmpi.back()].setPs2Energy( 0 );
00961 (*pfCandidates_)[tmpi.back()].addElementInBlock( blockref, kTrack.back() );
00962 }
00963 }
00964
00965
00966 if ( debug_ ) std::cout << "Loop over all associated ECAL clusters" << std::endl;
00967
00968 for(IE ie = ecalElems.begin(); ie != ecalElems.end(); ++ie ) {
00969
00970 unsigned index = ie->second;
00971 PFBlockElement::Type type = elements[index].type();
00972 assert( type == PFBlockElement::ECAL );
00973 if ( debug_ ) std::cout << elements[index] << std::endl;
00974
00975
00976 if ( debug_ && ! active[index] ) std::cout << "is not active - ignore " << std::endl;
00977 if ( ! active[index] ) continue;
00978
00979
00980
00981 block.associatedElements( index, linkData,
00982 sortedTracks,
00983 reco::PFBlockElement::TRACK,
00984 reco::PFBlock::LINKTEST_ALL );
00985 bool skip = true;
00986 for (unsigned ic=0; ic<kTrack.size();++ic) {
00987 if ( sortedTracks.begin()->second == kTrack[ic] ) {
00988 skip = false;
00989 break;
00990 }
00991 }
00992 if ( debug_ && skip ) std::cout << "is closer to another track - ignore " << std::endl;
00993 if ( skip ) continue;
00994
00995
00996
00997 reco::PFClusterRef clusterRef = elements[index].clusterRef();
00998 assert( !clusterRef.isNull() );
00999
01000 if ( debug_ ) {
01001 double dist = ie->first;
01002 std::cout <<"Ecal cluster with raw energy = " << clusterRef->energy()
01003 <<" linked with distance = " << dist << std::endl;
01004 }
01005
01006
01007
01008
01009
01010
01011
01012
01013
01014
01015
01016 vector<double> ps1Ene(1,static_cast<double>(0.));
01017 associatePSClusters(index, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
01018 vector<double> ps2Ene(1,static_cast<double>(0.));
01019 associatePSClusters(index, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
01020
01021
01022 bool crackCorrection = false;
01023 double ecalEnergy = calibration_->energyEm(*clusterRef,ps1Ene,ps2Ene,crackCorrection);
01024 if ( debug_ )
01025 std::cout << "Corrected ECAL(+PS) energy = " << ecalEnergy << std::endl;
01026
01027
01028
01029 totalEcal += ecalEnergy;
01030 totalHcal = 0.;
01031 double previousCalibEcal = calibEcal;
01032 double previousSlopeEcal = slopeEcal;
01033 calibEcal = std::max(totalEcal,0.);
01034 calibHcal = 0.;
01035 calibration_->energyEmHad(trackMomentum,calibEcal,calibHcal,
01036 clusterRef->positionREP().Eta(),
01037 clusterRef->positionREP().Phi());
01038 if ( totalEcal > 0.) slopeEcal = calibEcal/totalEcal;
01039
01040 if ( debug_ )
01041 std::cout << "The total calibrated energy so far amounts to = " << calibEcal << std::endl;
01042
01043
01044
01045 if ( connectedToEcal && calibEcal - trackMomentum >= 0. ) {
01046
01047
01048 calibEcal = previousCalibEcal;
01049 slopeEcal = previousSlopeEcal;
01050 totalEcal = calibEcal/slopeEcal;
01051
01052
01053
01054 active[index] = false;
01055
01056
01057 std::multimap<double, unsigned> assTracks;
01058 block.associatedElements( index, linkData,
01059 assTracks,
01060 reco::PFBlockElement::TRACK,
01061 reco::PFBlock::LINKTEST_ALL );
01062
01063
01064 unsigned tmpe = reconstructCluster( *clusterRef, ecalEnergy );
01065 (*pfCandidates_)[tmpe].setEcalEnergy( clusterRef->energy(), ecalEnergy );
01066 (*pfCandidates_)[tmpe].setHcalEnergy( 0., 0. );
01067 (*pfCandidates_)[tmpe].setPs1Energy( ps1Ene[0] );
01068 (*pfCandidates_)[tmpe].setPs2Energy( ps2Ene[0] );
01069 (*pfCandidates_)[tmpe].addElementInBlock( blockref, index );
01070 (*pfCandidates_)[tmpe].addElementInBlock( blockref, assTracks.begin()->second );
01071 break;
01072 }
01073
01074
01075 connectedToEcal = true;
01076 iEcal = index;
01077 active[index] = false;
01078 for (unsigned ic=0; ic<tmpi.size();++ic)
01079 (*pfCandidates_)[tmpi[ic]].addElementInBlock( blockref, iEcal );
01080
01081
01082 }
01083
01084 bool bNeutralProduced = false;
01085
01086
01087 if( connectedToEcal ) {
01088 reco::PFClusterRef pivotalRef = elements[iEcal].clusterRef();
01089
01090 double neutralEnergy = calibEcal - trackMomentum;
01091
01092
01093
01094
01095
01096
01097
01098
01099
01100
01101
01102
01103
01104
01105
01106
01107
01108
01109
01110
01111
01112
01113
01114
01115
01116
01117
01118
01119
01120
01121
01122
01123
01124
01125
01126 double resol = neutralHadronEnergyResolution(trackMomentum,pivotalRef->positionREP().Eta());
01127 resol *= trackMomentum;
01128 if ( neutralEnergy > std::max(0.5,nSigmaECAL_*resol) ) {
01129 neutralEnergy /= slopeEcal;
01130 unsigned tmpj = reconstructCluster( *pivotalRef, neutralEnergy );
01131 (*pfCandidates_)[tmpj].setEcalEnergy( pivotalRef->energy(), neutralEnergy );
01132 (*pfCandidates_)[tmpj].setHcalEnergy( 0., 0. );
01133 (*pfCandidates_)[tmpj].setPs1Energy( 0. );
01134 (*pfCandidates_)[tmpj].setPs2Energy( 0. );
01135 (*pfCandidates_)[tmpj].addElementInBlock(blockref, iEcal);
01136 bNeutralProduced = true;
01137 for (unsigned ic=0; ic<kTrack.size();++ic)
01138 (*pfCandidates_)[tmpj].addElementInBlock( blockref, kTrack[ic] );
01139 }
01140
01141
01142 for (unsigned ic=0; ic<tmpi.size();++ic) {
01143
01144
01145 if ( (*pfCandidates_)[tmpi[ic]].particleId() == reco::PFCandidate::mu ) continue;
01146
01147 double fraction = (*pfCandidates_)[tmpi[ic]].trackRef()->p()/trackMomentum;
01148 double ecalCal = bNeutralProduced ?
01149 (calibEcal-neutralEnergy*slopeEcal)*fraction : calibEcal*fraction;
01150 double ecalRaw = totalEcal*fraction;
01151
01152 if (debug_) cout << "The fraction after photon supression is " << fraction << " calibrated ecal = " << ecalCal << endl;
01153
01154 (*pfCandidates_)[tmpi[ic]].setEcalEnergy( ecalRaw, ecalCal );
01155 (*pfCandidates_)[tmpi[ic]].setHcalEnergy( 0., 0. );
01156 (*pfCandidates_)[tmpi[ic]].setPs1Energy( 0 );
01157 (*pfCandidates_)[tmpi[ic]].setPs2Energy( 0 );
01158 (*pfCandidates_)[tmpi[ic]].addElementInBlock( blockref, kTrack[ic] );
01159 }
01160
01161 }
01162
01163
01164 for (unsigned ic=0; ic<tmpi.size();++ic) {
01165 const PFCandidate& pfc = (*pfCandidates_)[tmpi[ic]];
01166 const PFCandidate::ElementsInBlocks& eleInBlocks = pfc.elementsInBlocks();
01167 if ( eleInBlocks.size() == 0 ) {
01168 if ( debug_ )std::cout << "Single track / Fill element in block! " << std::endl;
01169 (*pfCandidates_)[tmpi[ic]].addElementInBlock( blockref, kTrack[ic] );
01170 }
01171 }
01172
01173 }
01174
01175
01176
01177
01178 for(IE ie = hcalElems.begin(); ie != hcalElems.end(); ++ie ) {
01179
01180 unsigned index = ie->second;
01181
01182 PFBlockElement::Type type = elements[index].type();
01183
01184 if(debug_) {
01185 double dist = block.dist(iTrack,index,linkData,reco::PFBlock::LINKTEST_ALL);
01186 cout<<"\telement "<<elements[index]<<" linked with distance "<< dist <<endl;
01187 }
01188
01189 assert( type == PFBlockElement::HCAL );
01190
01191
01192
01193 if( !hcalFound ) {
01194 if(debug_)
01195 cout<<"\t\tclosest hcal cluster, doing nothing"<<endl;
01196
01197 hcalFound = true;
01198
01199
01200
01201 }
01202 else {
01203
01204 if(debug_)
01205 cout<<"\t\tsecondary hcal cluster. unlinking"<<endl;
01206 block.setLink( iTrack, index, -1., linkData,
01207 PFBlock::LINKTEST_RECHIT );
01208 }
01209 }
01210 }
01211
01212
01213
01214
01215 if( !(hfEmIs.empty() && hfHadIs.empty() ) ) {
01216
01217
01218 assert( hfEmIs.size() + hfHadIs.size() == elements.size() );
01219
01220 if( elements.size() == 1 ) {
01221
01222 reco::PFClusterRef clusterRef = elements[0].clusterRef();
01223 double energyHF = 0.;
01224 double uncalibratedenergyHF = 0.;
01225 unsigned tmpi = 0;
01226 switch( clusterRef->layer() ) {
01227 case PFLayer::HF_EM:
01228
01229 energyHF = clusterRef->energy();
01230 uncalibratedenergyHF = energyHF;
01231 if(thepfEnergyCalibrationHF_->getcalibHF_use()==true){
01232 energyHF = thepfEnergyCalibrationHF_->energyEm(uncalibratedenergyHF,
01233 clusterRef->positionREP().Eta(),
01234 clusterRef->positionREP().Phi());
01235 }
01236 tmpi = reconstructCluster( *clusterRef, energyHF );
01237 (*pfCandidates_)[tmpi].setEcalEnergy( uncalibratedenergyHF, energyHF );
01238 (*pfCandidates_)[tmpi].setHcalEnergy( 0., 0.);
01239 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
01240 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
01241 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfEmIs[0] );
01242
01243 break;
01244 case PFLayer::HF_HAD:
01245
01246 energyHF = clusterRef->energy();
01247 uncalibratedenergyHF = energyHF;
01248 if(thepfEnergyCalibrationHF_->getcalibHF_use()==true){
01249 energyHF = thepfEnergyCalibrationHF_->energyHad(uncalibratedenergyHF,
01250 clusterRef->positionREP().Eta(),
01251 clusterRef->positionREP().Phi());
01252 }
01253 tmpi = reconstructCluster( *clusterRef, energyHF );
01254 (*pfCandidates_)[tmpi].setHcalEnergy( uncalibratedenergyHF, energyHF );
01255 (*pfCandidates_)[tmpi].setEcalEnergy( 0., 0.);
01256 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
01257 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
01258 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfHadIs[0] );
01259
01260 break;
01261 default:
01262 assert(0);
01263 }
01264 }
01265 else if( elements.size() == 2 ) {
01266
01267 reco::PFClusterRef c0 = elements[0].clusterRef();
01268 reco::PFClusterRef c1 = elements[1].clusterRef();
01269
01270 reco::PFClusterRef cem = c1;
01271 reco::PFClusterRef chad = c0;
01272
01273 if( c0->layer()== PFLayer::HF_EM ) {
01274 if ( c1->layer()== PFLayer::HF_HAD ) {
01275 cem = c0; chad = c1;
01276 } else {
01277 assert(0);
01278 }
01279
01280 double energyHfEm = cem->energy();
01281 double energyHfHad = chad->energy();
01282 double uncalibratedenergyHFEm = energyHfEm;
01283 double uncalibratedenergyHFHad = energyHfHad;
01284 if(thepfEnergyCalibrationHF_->getcalibHF_use()==true){
01285
01286 energyHfEm = thepfEnergyCalibrationHF_->energyEmHad(uncalibratedenergyHFEm,
01287 0.0,
01288 c0->positionREP().Eta(),
01289 c0->positionREP().Phi());
01290 energyHfHad = thepfEnergyCalibrationHF_->energyEmHad(0.0,
01291 uncalibratedenergyHFHad,
01292 c1->positionREP().Eta(),
01293 c1->positionREP().Phi());
01294 }
01295 unsigned tmpi = reconstructCluster( *chad, energyHfEm+energyHfHad );
01296 (*pfCandidates_)[tmpi].setEcalEnergy( uncalibratedenergyHFEm, energyHfEm );
01297 (*pfCandidates_)[tmpi].setHcalEnergy( uncalibratedenergyHFHad, energyHfHad);
01298 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
01299 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
01300 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfEmIs[0] );
01301 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfHadIs[0] );
01302
01303
01304 }
01305 else {
01306 cerr<<"Warning: 2 elements, but not 1 HFEM and 1 HFHAD"<<endl;
01307 cerr<<block<<endl;
01308
01309
01310 }
01311 }
01312 else {
01313
01314
01315 cerr<<"Warning: HF, but n elem different from 1 or 2"<<endl;
01316 cerr<<block<<endl;
01317
01318
01319 }
01320 }
01321
01322
01323
01324 if(debug_) {
01325 cout<<endl;
01326 cout<<endl<<"--------------- loop hcal ---------------------"<<endl;
01327 }
01328
01329
01330
01331
01332
01333
01334
01335 for(unsigned i=0; i<hcalIs.size(); i++) {
01336
01337 unsigned iHcal= hcalIs[i];
01338 PFBlockElement::Type type = elements[iHcal].type();
01339
01340 assert( type == PFBlockElement::HCAL );
01341
01342 if(debug_) cout<<endl<<elements[iHcal]<<endl;
01343
01344
01345
01346
01347
01348 std::multimap<double, unsigned> sortedTracks;
01349 block.associatedElements( iHcal, linkData,
01350 sortedTracks,
01351 reco::PFBlockElement::TRACK,
01352 reco::PFBlock::LINKTEST_ALL );
01353
01354 std::multimap< unsigned, std::pair<double, unsigned> > associatedEcals;
01355
01356 std::map< unsigned, std::pair<double, double> > associatedPSs;
01357
01358 std::multimap<double, std::pair<unsigned,bool> > associatedTracks;
01359
01360
01361 std::multimap<double,std::pair<unsigned,math::XYZVector> > ecalSatellites;
01362 std::pair<unsigned,math::XYZVector> fakeSatellite = make_pair(iHcal,math::XYZVector(0.,0.,0.));
01363 ecalSatellites.insert( make_pair(-1., fakeSatellite) );
01364
01365 PFClusterRef hclusterref = elements[iHcal].clusterRef();
01366 assert(!hclusterref.isNull() );
01367
01368
01369
01370
01371
01372
01373
01374
01375
01376 if( sortedTracks.empty() ) {
01377 if(debug_)
01378 cout<<"\tno associated tracks, keep for later"<<endl;
01379 continue;
01380 }
01381
01382
01383 active[iHcal] = false;
01384
01385
01386
01387
01388
01389
01390
01391 if(debug_) cout<<"\t"<<sortedTracks.size()<<" associated tracks:"<<endl;
01392
01393 double totalChargedMomentum = 0;
01394 double sumpError2 = 0.;
01395 double totalEcal = 0.;
01396 double totalHcal = hclusterref->energy();
01397 vector<double> hcalP;
01398 vector<double> hcalDP;
01399 vector<unsigned> tkIs;
01400 double maxDPovP = -9999.;
01401
01402
01403 vector< unsigned > chargedHadronsIndices;
01404 vector< unsigned > chargedHadronsInBlock;
01405 double mergedNeutralHadronEnergy = 0;
01406 double mergedPhotonEnergy = 0;
01407 double muonHCALEnergy = 0.;
01408 double muonECALEnergy = 0.;
01409 double muonHCALError = 0.;
01410 double muonECALError = 0.;
01411 unsigned nMuons = 0;
01412
01413
01414 math::XYZVector photonAtECAL(0.,0.,0.);
01415 math::XYZVector hadronDirection(hclusterref->position().X(),
01416 hclusterref->position().Y(),
01417 hclusterref->position().Z());
01418 hadronDirection = hadronDirection.Unit();
01419 math::XYZVector hadronAtECAL = totalHcal * hadronDirection;
01420
01421
01422 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
01423
01424 unsigned iTrack = ie->second;
01425
01426
01427 if ( !active[iTrack] ) continue;
01428
01429 PFBlockElement::Type type = elements[iTrack].type();
01430 assert( type == reco::PFBlockElement::TRACK );
01431
01432 reco::TrackRef trackRef = elements[iTrack].trackRef();
01433 assert( !trackRef.isNull() );
01434
01435
01436 const math::XYZPointF& chargedPosition =
01437 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[iTrack])->positionAtECALEntrance();
01438 math::XYZVector chargedDirection(chargedPosition.X(),chargedPosition.Y(),chargedPosition.Z());
01439 chargedDirection = chargedDirection.Unit();
01440
01441
01442 std::multimap<double, unsigned> sortedEcals;
01443 block.associatedElements( iTrack, linkData,
01444 sortedEcals,
01445 reco::PFBlockElement::ECAL,
01446 reco::PFBlock::LINKTEST_ALL );
01447
01448 if(debug_) cout<<"\t\t\tnumber of Ecal elements linked to this track: "
01449 <<sortedEcals.size()<<endl;
01450
01451
01452 reco::MuonRef muonRef = elements[iTrack].muonRef();
01453 bool thisIsAMuon = PFMuonAlgo::isMuon(elements[iTrack]);
01454 bool thisIsAnIsolatedMuon = PFMuonAlgo::isIsolatedMuon(elements[iTrack]);
01455 bool thisIsALooseMuon = false;
01456
01457 if(!thisIsAMuon ) {
01458 thisIsALooseMuon = PFMuonAlgo::isLooseMuon(elements[iTrack]);
01459 }
01460 if ( thisIsAMuon ) {
01461 if ( debug_ ) {
01462 std::cout << "\t\tThis track is identified as a muon - remove it from the stack" << std::endl;
01463 std::cout << "\t\t" << elements[iTrack] << std::endl;
01464 }
01465
01466
01467 unsigned tmpi = reconstructTrack( elements[iTrack] );
01468 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
01469 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
01470 double muonHcal = std::min(muonHCAL_[0]+muonHCAL_[1],totalHcal);
01471
01472
01473 bool letMuonEatCaloEnergy = false;
01474
01475 if(thisIsAnIsolatedMuon){
01476
01477 double totalCaloEnergy = totalHcal / 1.30;
01478 unsigned iEcal = 0;
01479 if( !sortedEcals.empty() ) {
01480 iEcal = sortedEcals.begin()->second;
01481 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01482 totalCaloEnergy += eclusterref->energy();
01483 }
01484
01485
01486 if( (pfCandidates_->back()).p() > totalCaloEnergy ) letMuonEatCaloEnergy = true;
01487 }
01488
01489 if(letMuonEatCaloEnergy) muonHcal = totalHcal;
01490 double muonEcal =0.;
01491 unsigned iEcal = 0;
01492 if( !sortedEcals.empty() ) {
01493 iEcal = sortedEcals.begin()->second;
01494 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01495 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal);
01496 muonEcal = std::min(muonECAL_[0]+muonECAL_[1],eclusterref->energy());
01497 if(letMuonEatCaloEnergy) muonEcal = eclusterref->energy();
01498
01499 if ( eclusterref->energy() - muonEcal < 0.2 ) active[iEcal] = false;
01500 (*pfCandidates_)[tmpi].setEcalEnergy(eclusterref->energy(), muonEcal);
01501 }
01502 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal,muonHcal);
01503
01504 if(letMuonEatCaloEnergy){
01505 muonHCALEnergy += totalHcal;
01506 muonHCALError += 0.;
01507 muonECALEnergy += muonEcal;
01508 muonECALError += 0.;
01509 photonAtECAL -= muonEcal*chargedDirection;
01510 hadronAtECAL -= totalHcal*chargedDirection;
01511 active[iEcal] = false;
01512 active[iHcal] = false; }
01513 else{
01514
01515 muonHCALEnergy += muonHCAL_[0];
01516 muonHCALError += muonHCAL_[1]*muonHCAL_[1];
01517 muonECALEnergy += muonECAL_[0];
01518 muonECALError += muonECAL_[1]*muonECAL_[1];
01519
01520 photonAtECAL -= muonECAL_[0]*chargedDirection;
01521 hadronAtECAL -= muonHCAL_[0]*chargedDirection;
01522 }
01523
01524
01525 active[iTrack] = false;
01526
01527 continue;
01528 }
01529
01530
01531
01532 if(debug_) cout<<"\t\t"<<elements[iTrack]<<endl;
01533
01534
01535 double trackMomentum = trackRef->p();
01536 totalChargedMomentum += trackMomentum;
01537
01538
01539
01540 if ( thisIsALooseMuon && !thisIsAMuon ) nMuons += 1;
01541
01542
01543
01544
01545 double Dpt = trackRef->ptError();
01546 double blowError = 1.;
01547 switch (trackRef->algo()) {
01548 case TrackBase::ctf:
01549 case TrackBase::iter0:
01550 case TrackBase::iter1:
01551 case TrackBase::iter2:
01552 case TrackBase::iter3:
01553 blowError = 1.;
01554 break;
01555 case TrackBase::iter4:
01556 blowError = factors45_[0];
01557 break;
01558 case TrackBase::iter5:
01559 blowError = factors45_[1];
01560 break;
01561 default:
01562 blowError = 1E9;
01563 break;
01564 }
01565
01566 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
01567
01568 if ( isPrimaryOrSecondary ) blowError = 1.;
01569
01570 std::pair<unsigned,bool> tkmuon = make_pair(iTrack,thisIsALooseMuon);
01571 associatedTracks.insert( make_pair(-Dpt*blowError, tkmuon) );
01572
01573
01574 double Dp = trackRef->qoverpError()*trackMomentum*trackMomentum;
01575 sumpError2 += Dp*Dp;
01576
01577 bool connectedToEcal = false;
01578 if( !sortedEcals.empty() )
01579 {
01580
01581
01582 for ( IE iec=sortedEcals.begin();
01583 iec!=sortedEcals.end(); ++iec ) {
01584
01585 unsigned iEcal = iec->second;
01586 double dist = iec->first;
01587
01588
01589 if( !active[iEcal] ) {
01590 if(debug_) cout<<"\t\tcluster locked"<<endl;
01591 continue;
01592 }
01593
01594
01595 PFBlockElement::Type type = elements[ iEcal ].type();
01596 assert( type == PFBlockElement::ECAL );
01597 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01598 assert(!eclusterref.isNull() );
01599
01600
01601 std::multimap<double, unsigned> sortedTracksEcal;
01602 block.associatedElements( iEcal, linkData,
01603 sortedTracksEcal,
01604 reco::PFBlockElement::TRACK,
01605 reco::PFBlock::LINKTEST_ALL );
01606 unsigned jTrack = sortedTracksEcal.begin()->second;
01607 if ( jTrack != iTrack ) continue;
01608
01609
01610
01611 double distEcal = block.dist(jTrack,iEcal,linkData,
01612 reco::PFBlock::LINKTEST_ALL);
01613
01614
01615
01616
01617
01618
01619 vector<double> ps1Ene(1,static_cast<double>(0.));
01620 associatePSClusters(iEcal, reco::PFBlockElement::PS1,
01621 block, elements, linkData, active,
01622 ps1Ene);
01623 vector<double> ps2Ene(1,static_cast<double>(0.));
01624 associatePSClusters(iEcal, reco::PFBlockElement::PS2,
01625 block, elements, linkData, active,
01626 ps2Ene);
01627 std::pair<double,double> psEnes = make_pair(ps1Ene[0],
01628 ps2Ene[0]);
01629 associatedPSs.insert(make_pair(iEcal,psEnes));
01630
01631
01632 bool crackCorrection = false;
01633 float ecalEnergyCalibrated = calibration_->energyEm(*eclusterref,ps1Ene,ps2Ene,crackCorrection);
01634 math::XYZVector photonDirection(eclusterref->position().X(),
01635 eclusterref->position().Y(),
01636 eclusterref->position().Z());
01637 photonDirection = photonDirection.Unit();
01638
01639 if ( !connectedToEcal ) {
01640
01641 if(debug_) cout<<"\t\t\tclosest: "
01642 <<elements[iEcal]<<endl;
01643
01644 connectedToEcal = true;
01645
01646
01647
01648 std::pair<unsigned,math::XYZVector> satellite =
01649 make_pair(iEcal,ecalEnergyCalibrated*photonDirection);
01650 ecalSatellites.insert( make_pair(-1., satellite) );
01651
01652 } else {
01653
01654 std::pair<unsigned,math::XYZVector> satellite =
01655 make_pair(iEcal,ecalEnergyCalibrated*photonDirection);
01656 ecalSatellites.insert( make_pair(dist, satellite) );
01657
01658 }
01659
01660 std::pair<double, unsigned> associatedEcal
01661 = make_pair( distEcal, iEcal );
01662 associatedEcals.insert( make_pair(iTrack, associatedEcal) );
01663
01664 }
01665 }
01666
01667 }
01668
01669
01670
01671 double caloEnergy = 0.;
01672 double slopeEcal = 1.0;
01673 double calibEcal = 0.;
01674 double calibHcal = 0.;
01675 hadronDirection = hadronAtECAL.Unit();
01676
01677
01678 double Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
01679 Caloresolution *= totalChargedMomentum;
01680
01681 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
01682 totalEcal -= std::min(totalEcal,muonECALEnergy);
01683 totalHcal -= std::min(totalHcal,muonHCALEnergy);
01684 if ( totalEcal < 1E-9 ) photonAtECAL = math::XYZVector(0.,0.,0.);
01685 if ( totalHcal < 1E-9 ) hadronAtECAL = math::XYZVector(0.,0.,0.);
01686
01687
01688
01689
01690
01691 for ( IS is = ecalSatellites.begin(); is != ecalSatellites.end(); ++is ) {
01692
01693
01694 double previousCalibEcal = calibEcal;
01695 double previousCalibHcal = calibHcal;
01696 double previousCaloEnergy = caloEnergy;
01697 double previousSlopeEcal = slopeEcal;
01698 math::XYZVector previousHadronAtECAL = hadronAtECAL;
01699
01700 totalEcal += sqrt(is->second.second.Mag2());
01701 photonAtECAL += is->second.second;
01702 calibEcal = std::max(0.,totalEcal);
01703 calibHcal = std::max(0.,totalHcal);
01704 hadronAtECAL = calibHcal * hadronDirection;
01705
01706
01707 calibration_->energyEmHad(totalChargedMomentum,calibEcal,calibHcal,
01708 hclusterref->positionREP().Eta(),
01709 hclusterref->positionREP().Phi());
01710 caloEnergy = calibEcal+calibHcal;
01711 if ( totalEcal > 0.) slopeEcal = calibEcal/totalEcal;
01712
01713 hadronAtECAL = calibHcal * hadronDirection;
01714
01715
01716
01717 if ( is->first < 0. || caloEnergy - totalChargedMomentum <= 0. ) {
01718 if(debug_) cout<<"\t\t\tactive, adding "<<is->second.second
01719 <<" to ECAL energy, and locking"<<endl;
01720 active[is->second.first] = false;
01721 continue;
01722 }
01723
01724
01725
01726 totalEcal -= sqrt(is->second.second.Mag2());
01727 photonAtECAL -= is->second.second;
01728 calibEcal = previousCalibEcal;
01729 calibHcal = previousCalibHcal;
01730 hadronAtECAL = previousHadronAtECAL;
01731 slopeEcal = previousSlopeEcal;
01732 caloEnergy = previousCaloEnergy;
01733
01734 break;
01735
01736 }
01737
01738
01739 assert(caloEnergy>=0);
01740
01741
01742
01743
01744
01745
01746
01747
01748
01749
01750
01751
01752 double TotalError = sqrt(sumpError2 + Caloresolution*Caloresolution);
01753
01754
01756 if ( totalChargedMomentum - caloEnergy > nSigmaTRACK_*Caloresolution ) {
01757
01758
01759
01760
01761
01762
01763
01764
01765
01766
01767
01768
01769
01770 if ( nMuons > 0 ) {
01771 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
01772 unsigned iTrack = it->second.first;
01773
01774 if ( !active[iTrack] ) continue;
01775
01776 if ( !it->second.second ) continue;
01777
01778 bool isTrack = elements[it->second.first].muonRef()->isTrackerMuon();
01779 double trackMomentum = elements[it->second.first].trackRef()->p();
01780
01781
01782 std::multimap<double, unsigned> sortedEcals;
01783 block.associatedElements( iTrack, linkData,
01784 sortedEcals,
01785 reco::PFBlockElement::ECAL,
01786 reco::PFBlock::LINKTEST_ALL );
01787
01788 unsigned tmpi = reconstructTrack( elements[iTrack] );
01789 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
01790 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
01791 double muonHcal = std::min(muonHCAL_[0]+muonHCAL_[1],totalHcal);
01792 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal,muonHcal);
01793 if( !sortedEcals.empty() ) {
01794 unsigned iEcal = sortedEcals.begin()->second;
01795 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01796 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal);
01797 double muonEcal = std::min(muonECAL_[0]+muonECAL_[1],eclusterref->energy());
01798 (*pfCandidates_)[tmpi].setEcalEnergy(eclusterref->energy(),muonEcal);
01799 }
01800 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::mu;
01801 (*pfCandidates_)[tmpi].setParticleType(particleType);
01802
01803 if(PFMuonAlgo::isGlobalLooseMuon(elements[it->second.first].muonRef())){
01804
01805 double muonMomentum = elements[it->second.first].muonRef()->combinedMuon()->p();
01806 double staMomentum = elements[it->second.first].muonRef()->standAloneMuon()->p();
01807 double trackPtError = elements[it->second.first].trackRef()->ptError();
01808 double muonPtError = elements[it->second.first].muonRef()->combinedMuon()->ptError();
01809 double staPtError = elements[it->second.first].muonRef()->standAloneMuon()->ptError();
01810 double globalCorr = 1;
01811
01812
01813 if(!isTrack || muonPtError < trackPtError || staPtError < trackPtError){
01814 if(muonPtError < staPtError) globalCorr = muonMomentum/trackMomentum;
01815 else globalCorr = staMomentum/trackMomentum;
01816 }
01817 (*pfCandidates_)[tmpi].rescaleMomentum(globalCorr);
01818
01819 if (debug_){
01820 std::cout << "\tElement " << elements[iTrack] << std::endl
01821 << "PFAlgo: particle type set to muon (global, loose), pT = "<<elements[it->second.first].muonRef()->pt()<<std::endl;
01822 PFMuonAlgo::printMuonProperties(elements[it->second.first].muonRef());
01823 }
01824 }
01825 else{
01826 if (debug_){
01827 std::cout << "\tElement " << elements[iTrack] << std::endl
01828 << "PFAlgo: particle type set to muon (tracker, loose), pT = "<<elements[it->second.first].muonRef()->pt()<<std::endl;
01829 PFMuonAlgo::printMuonProperties(elements[it->second.first].muonRef());
01830 }
01831 }
01832
01833
01834 const math::XYZPointF& chargedPosition =
01835 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[it->second.first])->positionAtECALEntrance();
01836 math::XYZVector chargedDirection(chargedPosition.X(), chargedPosition.Y(), chargedPosition.Z());
01837 chargedDirection = chargedDirection.Unit();
01838 totalChargedMomentum -= trackMomentum;
01839
01840 if ( totalEcal > 0. )
01841 calibEcal -= std::min(calibEcal,muonECAL_[0]*calibEcal/totalEcal);
01842 if ( totalHcal > 0. )
01843 calibHcal -= std::min(calibHcal,muonHCAL_[0]*calibHcal/totalHcal);
01844 totalEcal -= std::min(totalEcal,muonECAL_[0]);
01845 totalHcal -= std::min(totalHcal,muonHCAL_[0]);
01846 if ( totalEcal > muonECAL_[0] ) photonAtECAL -= muonECAL_[0] * chargedDirection;
01847 if ( totalHcal > muonHCAL_[0] ) hadronAtECAL -= muonHCAL_[0]*calibHcal/totalHcal * chargedDirection;
01848 caloEnergy = calibEcal+calibHcal;
01849 muonHCALEnergy += muonHCAL_[0];
01850 muonHCALError += muonHCAL_[1]*muonHCAL_[1];
01851 muonECALEnergy += muonECAL_[0];
01852 muonECALError += muonECAL_[1]*muonECAL_[1];
01853 active[iTrack] = false;
01854
01855
01856
01857 }
01858
01859 Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
01860 Caloresolution *= totalChargedMomentum;
01861 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
01862 }
01863 }
01864
01865
01866
01867
01868
01869
01870
01871
01872
01873
01874
01875
01876
01877
01878
01879 unsigned corrTrack = 10000000;
01880 double corrFact = 1.;
01881
01882 if (rejectTracks_Bad_ &&
01883 totalChargedMomentum - caloEnergy > nSigmaTRACK_*Caloresolution) {
01884
01885 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
01886 unsigned iTrack = it->second.first;
01887
01888 if ( !active[iTrack] ) continue;
01889 const reco::TrackRef& trackref = elements[it->second.first].trackRef();
01890
01891 double dptRel = fabs(it->first)/trackref->pt()*100;
01892 bool isSecondary = isFromSecInt(elements[iTrack], "secondary");
01893 bool isPrimary = isFromSecInt(elements[iTrack], "primary");
01894
01895 if ( isSecondary && dptRel < dptRel_DispVtx_ ) continue;
01896
01897 if ( fabs(it->first) < ptError_ ) break;
01898
01899 double wouldBeTotalChargedMomentum =
01900 totalChargedMomentum - trackref->p();
01901
01902
01903
01904 if ( wouldBeTotalChargedMomentum > caloEnergy ) {
01905
01906 if (debug_ && isSecondary) {
01907 cout << "In bad track rejection step dptRel = " << dptRel << " dptRel_DispVtx_ = " << dptRel_DispVtx_ << endl;
01908 cout << "The calo energy would be still smaller even without this track but it is attached to a NI"<< endl;
01909 }
01910
01911
01912 if(isPrimary || (isSecondary && dptRel < dptRel_DispVtx_)) continue;
01913 active[iTrack] = false;
01914 totalChargedMomentum = wouldBeTotalChargedMomentum;
01915 if ( debug_ )
01916 std::cout << "\tElement " << elements[iTrack]
01917 << " rejected (Dpt = " << -it->first
01918 << " GeV/c, algo = " << trackref->algo() << ")" << std::endl;
01919
01920 } else {
01921 if(isPrimary) break;
01922 corrTrack = iTrack;
01923 corrFact = (caloEnergy - wouldBeTotalChargedMomentum)/elements[it->second.first].trackRef()->p();
01924 if ( trackref->p()*corrFact < 0.05 ) {
01925 corrFact = 0.;
01926 active[iTrack] = false;
01927 }
01928 totalChargedMomentum -= trackref->p()*(1.-corrFact);
01929 if ( debug_ )
01930 std::cout << "\tElement " << elements[iTrack]
01931 << " (Dpt = " << -it->first
01932 << " GeV/c, algo = " << trackref->algo()
01933 << ") rescaled by " << corrFact
01934 << " Now the total charged momentum is " << totalChargedMomentum << endl;
01935 break;
01936 }
01937 }
01938 }
01939
01940
01941 Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
01942 Caloresolution *= totalChargedMomentum;
01943 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
01944
01945
01946
01947
01948 if ( rejectTracks_Step45_ &&
01949 sortedTracks.size() > 1 &&
01950 totalChargedMomentum - caloEnergy > nSigmaTRACK_*Caloresolution ) {
01951 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
01952 unsigned iTrack = it->second.first;
01953 reco::TrackRef trackref = elements[iTrack].trackRef();
01954 if ( !active[iTrack] ) continue;
01955
01956 double dptRel = fabs(it->first)/trackref->pt()*100;
01957 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
01958
01959
01960
01961
01962 if (isPrimaryOrSecondary && dptRel < dptRel_DispVtx_) continue;
01963
01964 switch (trackref->algo()) {
01965 case TrackBase::ctf:
01966 case TrackBase::iter0:
01967 case TrackBase::iter1:
01968 case TrackBase::iter2:
01969 case TrackBase::iter3:
01970 break;
01971 case TrackBase::iter4:
01972 case TrackBase::iter5:
01973 active[iTrack] = false;
01974 totalChargedMomentum -= trackref->p();
01975
01976 if ( debug_ )
01977 std::cout << "\tElement " << elements[iTrack]
01978 << " rejected (Dpt = " << -it->first
01979 << " GeV/c, algo = " << trackref->algo() << ")" << std::endl;
01980 break;
01981 default:
01982 break;
01983 }
01984 }
01985 }
01986
01987
01988 Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
01989 Caloresolution *= totalChargedMomentum;
01990 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
01991
01992
01993 sumpError2 = 0.;
01994 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
01995 unsigned iTrack = it->second.first;
01996 if ( !active[iTrack] ) continue;
01997 reco::TrackRef trackRef = elements[iTrack].trackRef();
01998 double trackMomentum = trackRef->p();
01999 double Dp = trackRef->qoverpError()*trackMomentum*trackMomentum;
02000 unsigned tmpi = reconstructTrack( elements[iTrack] );
02001 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
02002 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02003 std::pair<II,II> myEcals = associatedEcals.equal_range(iTrack);
02004 for (II ii=myEcals.first; ii!=myEcals.second; ++ii ) {
02005 unsigned iEcal = ii->second.second;
02006 if ( active[iEcal] ) continue;
02007 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02008 }
02009
02010 if ( iTrack == corrTrack ) {
02011 (*pfCandidates_)[tmpi].rescaleMomentum(corrFact);
02012 trackMomentum *= corrFact;
02013 }
02014 chargedHadronsIndices.push_back( tmpi );
02015 chargedHadronsInBlock.push_back( iTrack );
02016 active[iTrack] = false;
02017 hcalP.push_back(trackMomentum);
02018 hcalDP.push_back(Dp);
02019 if (Dp/trackMomentum > maxDPovP) maxDPovP = Dp/trackMomentum;
02020 sumpError2 += Dp*Dp;
02021 }
02022
02023
02024 TotalError = sqrt(sumpError2 + Caloresolution*Caloresolution);
02025
02026 if ( debug_ ) {
02027 cout<<"\tCompare Calo Energy to total charged momentum "<<endl;
02028 cout<<"\t\tsum p = "<<totalChargedMomentum<<" +- "<<sqrt(sumpError2)<<endl;
02029 cout<<"\t\tsum ecal = "<<totalEcal<<endl;
02030 cout<<"\t\tsum hcal = "<<totalHcal<<endl;
02031 cout<<"\t\t => Calo Energy = "<<caloEnergy<<" +- "<<Caloresolution<<endl;
02032 cout<<"\t\t => Calo Energy- total charged momentum = "
02033 <<caloEnergy-totalChargedMomentum<<" +- "<<TotalError<<endl;
02034 cout<<endl;
02035 }
02036
02037
02038
02040 double nsigma = nSigmaHCAL(totalChargedMomentum,hclusterref->positionREP().Eta());
02041
02042 if ( abs(totalChargedMomentum-caloEnergy)<nsigma*TotalError ) {
02043
02044
02045
02046
02047 if(debug_) {
02048 cout<<"\t\tcase 1: COMPATIBLE "
02049 <<"|Calo Energy- total charged momentum| = "
02050 <<abs(caloEnergy-totalChargedMomentum)
02051 <<" < "<<nsigma<<" x "<<TotalError<<endl;
02052 if (maxDPovP < 0.1 )
02053 cout<<"\t\t\tmax DP/P = "<< maxDPovP
02054 <<" less than 0.1: do nothing "<<endl;
02055 else
02056 cout<<"\t\t\tmax DP/P = "<< maxDPovP
02057 <<" > 0.1: take weighted averages "<<endl;
02058 }
02059
02060
02061 if (maxDPovP > 0.1) {
02062
02063
02064
02065 int nrows = chargedHadronsIndices.size();
02066 TMatrixTSym<double> a (nrows);
02067 TVectorD b (nrows);
02068 TVectorD check(nrows);
02069 double sigma2E = Caloresolution*Caloresolution;
02070 for(int i=0; i<nrows; i++) {
02071 double sigma2i = hcalDP[i]*hcalDP[i];
02072 if (debug_){
02073 cout<<"\t\t\ttrack associated to hcal "<<i
02074 <<" P = "<<hcalP[i]<<" +- "
02075 <<hcalDP[i]<<endl;
02076 }
02077 a(i,i) = 1./sigma2i + 1./sigma2E;
02078 b(i) = hcalP[i]/sigma2i + caloEnergy/sigma2E;
02079 for(int j=0; j<nrows; j++) {
02080 if (i==j) continue;
02081 a(i,j) = 1./sigma2E;
02082 }
02083 }
02084
02085
02086
02087
02088
02089
02090
02091
02092 TDecompChol decomp(a);
02093 bool ok = false;
02094 TVectorD x = decomp.Solve(b, ok);
02095
02096
02097 if (ok) {
02098 for (int i=0; i<nrows; i++){
02099
02100 unsigned ich = chargedHadronsIndices[i];
02101 double rescaleFactor = x(i)/hcalP[i];
02102 (*pfCandidates_)[ich].rescaleMomentum( rescaleFactor );
02103
02104 if(debug_){
02105 cout<<"\t\t\told p "<<hcalP[i]
02106 <<" new p "<<x(i)
02107 <<" rescale "<<rescaleFactor<<endl;
02108 }
02109 }
02110 }
02111 else {
02112 cerr<<"not ok!"<<endl;
02113 assert(0);
02114 }
02115 }
02116 }
02117
02119 else if( caloEnergy > totalChargedMomentum ) {
02120
02121
02122
02123
02124
02125
02126
02127
02128
02129
02130
02131
02132
02133
02134
02135
02136
02137
02138 double eNeutralHadron = caloEnergy - totalChargedMomentum;
02139 double ePhoton = (caloEnergy - totalChargedMomentum) / slopeEcal;
02140
02141 if(debug_) {
02142 if(!sortedTracks.empty() ){
02143 cout<<"\t\tcase 2: NEUTRAL DETECTION "
02144 <<caloEnergy<<" > "<<nsigma<<"x"<<TotalError
02145 <<" + "<<totalChargedMomentum<<endl;
02146 cout<<"\t\tneutral activity detected: "<<endl
02147 <<"\t\t\t photon = "<<ePhoton<<endl
02148 <<"\t\t\tor neutral hadron = "<<eNeutralHadron<<endl;
02149
02150 cout<<"\t\tphoton or hadron ?"<<endl;}
02151
02152 if(sortedTracks.empty() )
02153 cout<<"\t\tno track -> hadron "<<endl;
02154 else
02155 cout<<"\t\t"<<sortedTracks.size()
02156 <<"tracks -> check if the excess is photonic or hadronic"<<endl;
02157 }
02158
02159
02160 double ratioMax = 0.;
02161 reco::PFClusterRef maxEcalRef;
02162 unsigned maxiEcal= 9999;
02163
02164
02165
02166 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
02167
02168 unsigned iTrack = ie->second;
02169
02170 PFBlockElement::Type type = elements[iTrack].type();
02171 assert( type == reco::PFBlockElement::TRACK );
02172
02173 reco::TrackRef trackRef = elements[iTrack].trackRef();
02174 assert( !trackRef.isNull() );
02175
02176 II iae = associatedEcals.find(iTrack);
02177
02178 if( iae == associatedEcals.end() ) continue;
02179
02180
02181 unsigned iEcal = iae->second.second;
02182
02183 PFBlockElement::Type typeEcal = elements[iEcal].type();
02184 assert(typeEcal == reco::PFBlockElement::ECAL);
02185
02186 reco::PFClusterRef clusterRef = elements[iEcal].clusterRef();
02187 assert( !clusterRef.isNull() );
02188
02189 double pTrack = trackRef->p();
02190 double eECAL = clusterRef->energy();
02191 double eECALOverpTrack = eECAL / pTrack;
02192
02193 if ( eECALOverpTrack > ratioMax ) {
02194 ratioMax = eECALOverpTrack;
02195 maxEcalRef = clusterRef;
02196 maxiEcal = iEcal;
02197 }
02198
02199 }
02200
02201 std::vector<reco::PFClusterRef> pivotalClusterRef;
02202 std::vector<unsigned> iPivotal;
02203 std::vector<double> particleEnergy, ecalEnergy, hcalEnergy, rawecalEnergy, rawhcalEnergy;
02204 std::vector<math::XYZVector> particleDirection;
02205
02206
02207
02208 if ( ePhoton < totalEcal || eNeutralHadron-calibEcal < 1E-10 ) {
02209
02210 if ( !maxEcalRef.isNull() ) {
02211 particleEnergy.push_back(ePhoton);
02212 particleDirection.push_back(photonAtECAL);
02213 ecalEnergy.push_back(ePhoton);
02214 hcalEnergy.push_back(0.);
02215 rawecalEnergy.push_back(totalEcal);
02216 rawhcalEnergy.push_back(totalHcal);
02217 pivotalClusterRef.push_back(maxEcalRef);
02218 iPivotal.push_back(maxiEcal);
02219
02220 mergedPhotonEnergy = ePhoton;
02221 }
02222
02223 } else {
02224
02225
02226 if ( !maxEcalRef.isNull() ) {
02227 pivotalClusterRef.push_back(maxEcalRef);
02228 particleEnergy.push_back(totalEcal);
02229 particleDirection.push_back(photonAtECAL);
02230 ecalEnergy.push_back(totalEcal);
02231 hcalEnergy.push_back(0.);
02232 rawecalEnergy.push_back(totalEcal);
02233 rawhcalEnergy.push_back(totalHcal);
02234 iPivotal.push_back(maxiEcal);
02235
02236 mergedPhotonEnergy = totalEcal;
02237 }
02238
02239
02240 mergedNeutralHadronEnergy = eNeutralHadron-calibEcal;
02241 if ( mergedNeutralHadronEnergy > 1.0 ) {
02242 pivotalClusterRef.push_back(hclusterref);
02243 particleEnergy.push_back(mergedNeutralHadronEnergy);
02244 particleDirection.push_back(hadronAtECAL);
02245 ecalEnergy.push_back(0.);
02246 hcalEnergy.push_back(mergedNeutralHadronEnergy);
02247 rawecalEnergy.push_back(totalEcal);
02248 rawhcalEnergy.push_back(totalHcal);
02249 iPivotal.push_back(iHcal);
02250 }
02251
02252 }
02253
02254
02255
02256
02257
02258
02259 for ( unsigned iPivot=0; iPivot<iPivotal.size(); ++iPivot ) {
02260
02261 if ( particleEnergy[iPivot] < 0. )
02262 std::cout << "ALARM = Negative energy ! "
02263 << particleEnergy[iPivot] << std::endl;
02264
02265 bool useDirection = true;
02266 unsigned tmpi = reconstructCluster( *pivotalClusterRef[iPivot],
02267 particleEnergy[iPivot],
02268 useDirection,
02269 particleDirection[iPivot].X(),
02270 particleDirection[iPivot].Y(),
02271 particleDirection[iPivot].Z());
02272
02273
02274 (*pfCandidates_)[tmpi].setEcalEnergy( rawecalEnergy[iPivot],ecalEnergy[iPivot] );
02275 (*pfCandidates_)[tmpi].setHcalEnergy( rawhcalEnergy[iPivot],hcalEnergy[iPivot] );
02276 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
02277 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
02278 (*pfCandidates_)[tmpi].set_mva_nothing_gamma( -1. );
02279
02280
02281 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02282 for ( unsigned ich=0; ich<chargedHadronsInBlock.size(); ++ich) {
02283 unsigned iTrack = chargedHadronsInBlock[ich];
02284 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
02285 std::pair<II,II> myEcals = associatedEcals.equal_range(iTrack);
02286 for (II ii=myEcals.first; ii!=myEcals.second; ++ii ) {
02287 unsigned iEcal = ii->second.second;
02288 if ( active[iEcal] ) continue;
02289 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02290 }
02291 }
02292
02293 }
02294
02295 }
02296
02297
02298
02299
02300
02301 double totalHcalEnergyCalibrated = calibHcal;
02302 double totalEcalEnergyCalibrated = calibEcal;
02303
02304
02305
02306
02307
02308
02309
02310
02311
02312
02313
02314
02315
02316
02317 totalHcalEnergyCalibrated -= mergedNeutralHadronEnergy;
02318
02319 totalEcalEnergyCalibrated -= mergedPhotonEnergy;
02320
02321
02322
02323
02324 double chargedHadronsTotalEnergy = 0;
02325 for( unsigned ich=0; ich<chargedHadronsIndices.size(); ++ich ) {
02326 unsigned index = chargedHadronsIndices[ich];
02327 reco::PFCandidate& chargedHadron = (*pfCandidates_)[index];
02328 chargedHadronsTotalEnergy += chargedHadron.energy();
02329 }
02330
02331 for( unsigned ich=0; ich<chargedHadronsIndices.size(); ++ich ) {
02332 unsigned index = chargedHadronsIndices[ich];
02333 reco::PFCandidate& chargedHadron = (*pfCandidates_)[index];
02334 float fraction = chargedHadron.energy()/chargedHadronsTotalEnergy;
02335
02336 chargedHadron.setHcalEnergy( fraction * totalHcal, fraction * totalHcalEnergyCalibrated );
02337
02338 chargedHadron.setEcalEnergy( fraction * totalEcal, fraction * totalEcalEnergyCalibrated );
02339 }
02340
02341
02342 for ( IS is = ecalSatellites.begin(); is != ecalSatellites.end(); ++is ) {
02343
02344
02345 unsigned iEcal = is->second.first;
02346 if ( !active[iEcal] ) continue;
02347
02348
02349 PFBlockElement::Type type = elements[ iEcal ].type();
02350 assert( type == PFBlockElement::ECAL );
02351 PFClusterRef eclusterref = elements[iEcal].clusterRef();
02352 assert(!eclusterref.isNull() );
02353
02354
02355 active[iEcal] = false;
02356
02357
02358 std::multimap<double, unsigned> assTracks;
02359 block.associatedElements( iEcal, linkData,
02360 assTracks,
02361 reco::PFBlockElement::TRACK,
02362 reco::PFBlock::LINKTEST_ALL );
02363
02364
02365 unsigned tmpi = reconstructCluster( *eclusterref, sqrt(is->second.second.Mag2()) );
02366 (*pfCandidates_)[tmpi].setEcalEnergy( eclusterref->energy(),sqrt(is->second.second.Mag2()) );
02367 (*pfCandidates_)[tmpi].setHcalEnergy( 0., 0. );
02368 (*pfCandidates_)[tmpi].setPs1Energy( associatedPSs[iEcal].first );
02369 (*pfCandidates_)[tmpi].setPs2Energy( associatedPSs[iEcal].second );
02370 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02371 (*pfCandidates_)[tmpi].addElementInBlock( blockref, sortedTracks.begin()->second) ;
02372 }
02373
02374
02375 }
02376
02377
02378
02379
02380
02381 if(debug_) {
02382 cout<<endl;
02383 if(debug_)
02384 cout<<endl
02385 <<"---- loop remaining hcal ------- "<<endl;
02386 }
02387
02388
02389
02390
02391 for(unsigned ihcluster=0; ihcluster<hcalIs.size(); ihcluster++) {
02392 unsigned iHcal = hcalIs[ihcluster];
02393
02394 if(debug_)
02395 cout<<endl<<elements[iHcal]<<" ";
02396
02397
02398 if( !active[iHcal] ) {
02399 if(debug_)
02400 cout<<"not active"<<endl;
02401 continue;
02402 }
02403
02404
02405 std::multimap<double, unsigned> ecalElems;
02406 block.associatedElements( iHcal, linkData,
02407 ecalElems ,
02408 reco::PFBlockElement::ECAL,
02409 reco::PFBlock::LINKTEST_ALL );
02410
02411
02412 float totalEcal = 0.;
02413 float ecalMax = 0.;
02414 unsigned jEcal = 0;
02415 reco::PFClusterRef eClusterRef;
02416 for(IE ie = ecalElems.begin(); ie != ecalElems.end(); ++ie ) {
02417
02418 unsigned iEcal = ie->second;
02419 double dist = ie->first;
02420 PFBlockElement::Type type = elements[iEcal].type();
02421 assert( type == PFBlockElement::ECAL );
02422
02423
02424 if( !active[iEcal] ) continue;
02425
02426
02427 if ( dist > 0.15 ) continue;
02428
02429
02430
02431
02432
02433
02434
02435
02436
02437
02438
02439
02440
02441
02442
02443
02444
02445 std::multimap<double, unsigned> hcalElems;
02446 block.associatedElements( iEcal, linkData,
02447 hcalElems ,
02448 reco::PFBlockElement::HCAL,
02449 reco::PFBlock::LINKTEST_ALL );
02450
02451 bool isClosest = true;
02452 for(IE ih = hcalElems.begin(); ih != hcalElems.end(); ++ih ) {
02453
02454 unsigned jHcal = ih->second;
02455 double distH = ih->first;
02456
02457 if ( !active[jHcal] ) continue;
02458
02459 if ( distH < dist ) {
02460 isClosest = false;
02461 break;
02462 }
02463
02464 }
02465
02466 if (!isClosest) continue;
02467
02468
02469 if(debug_) {
02470 cout<<"\telement "<<elements[iEcal]<<" linked with dist "<< dist<<endl;
02471 cout<<"Added to HCAL cluster to form a neutral hadron"<<endl;
02472 }
02473
02474 reco::PFClusterRef eclusterRef = elements[iEcal].clusterRef();
02475 assert( !eclusterRef.isNull() );
02476
02477
02478 vector<double> ps1Ene(1,static_cast<double>(0.));
02479 associatePSClusters(iEcal, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
02480 vector<double> ps2Ene(1,static_cast<double>(0.));
02481 associatePSClusters(iEcal, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
02482 bool crackCorrection = false;
02483 double ecalEnergy = calibration_->energyEm(*eclusterRef,ps1Ene,ps2Ene,crackCorrection);
02484
02485
02486
02487 totalEcal += ecalEnergy;
02488 if ( ecalEnergy > ecalMax ) {
02489 ecalMax = ecalEnergy;
02490 eClusterRef = eclusterRef;
02491 jEcal = iEcal;
02492 }
02493
02494 active[iEcal] = false;
02495
02496 }
02497
02498
02499 PFClusterRef hclusterRef
02500 = elements[iHcal].clusterRef();
02501 assert( !hclusterRef.isNull() );
02502
02503
02504 double totalHcal = hclusterRef->energy();
02505
02506
02507
02508
02509 double caloEnergy = totalHcal;
02510
02511 double calibEcal = totalEcal > 0. ? totalEcal : 0.;
02512 double calibHcal = std::max(0.,totalHcal);
02513 if ( hclusterRef->layer() == PFLayer::HF_HAD ||
02514 hclusterRef->layer() == PFLayer::HF_EM ) {
02515 caloEnergy = totalHcal/0.7;
02516 calibEcal = totalEcal;
02517 calibHcal = caloEnergy;
02518 } else {
02519 calibration_->energyEmHad(-1.,calibEcal,calibHcal,
02520 hclusterRef->positionREP().Eta(),
02521 hclusterRef->positionREP().Phi());
02522 caloEnergy = calibEcal+calibHcal;
02523 }
02524
02525
02526
02527
02528
02529
02530
02531
02532
02533 unsigned tmpi = reconstructCluster( *hclusterRef,
02534 calibEcal+calibHcal );
02535
02536
02537 (*pfCandidates_)[tmpi].setEcalEnergy( totalEcal, calibEcal );
02538 (*pfCandidates_)[tmpi].setHcalEnergy( totalHcal, calibHcal );
02539 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
02540 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
02541 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02542
02543
02544 }
02545
02546
02547
02548
02549 if(debug_) {
02550 cout<<endl;
02551 if(debug_) cout<<endl<<"---- loop ecal------- "<<endl;
02552 }
02553
02554
02555
02556 for(unsigned i=0; i<ecalIs.size(); i++) {
02557 unsigned iEcal = ecalIs[i];
02558
02559 if(debug_)
02560 cout<<endl<<elements[iEcal]<<" ";
02561
02562 if( ! active[iEcal] ) {
02563 if(debug_)
02564 cout<<"not active"<<endl;
02565 continue;
02566 }
02567
02568 PFBlockElement::Type type = elements[ iEcal ].type();
02569 assert( type == PFBlockElement::ECAL );
02570
02571 PFClusterRef clusterref = elements[iEcal].clusterRef();
02572 assert(!clusterref.isNull() );
02573
02574 active[iEcal] = false;
02575
02576
02577
02578 vector<double> ps1Ene(1,static_cast<double>(0.));
02579 associatePSClusters(iEcal, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
02580 vector<double> ps2Ene(1,static_cast<double>(0.));
02581 associatePSClusters(iEcal, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
02582 bool crackCorrection = false;
02583 float ecalEnergy = calibration_->energyEm(*clusterref,ps1Ene,ps2Ene,crackCorrection);
02584
02585 double particleEnergy = ecalEnergy;
02586
02587 unsigned tmpi = reconstructCluster( *clusterref,
02588 particleEnergy );
02589
02590 (*pfCandidates_)[tmpi].setEcalEnergy( clusterref->energy(),ecalEnergy );
02591 (*pfCandidates_)[tmpi].setHcalEnergy( 0., 0. );
02592 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
02593 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
02594 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02595
02596
02597 }
02598
02599 }
02600
02601
02602
02603 unsigned PFAlgo::reconstructTrack( const reco::PFBlockElement& elt ) {
02604
02605 const reco::PFBlockElementTrack* eltTrack
02606 = dynamic_cast<const reco::PFBlockElementTrack*>(&elt);
02607
02608 reco::TrackRef trackRef = eltTrack->trackRef();
02609 const reco::Track& track = *trackRef;
02610
02611 reco::MuonRef muonRef = eltTrack->muonRef();
02612
02613 int charge = track.charge()>0 ? 1 : -1;
02614
02615
02616 double px = track.px();
02617 double py = track.py();
02618 double pz = track.pz();
02619 double energy = sqrt(track.p()*track.p() + 0.13957*0.13957);
02620
02621
02622
02623 bool thisIsAMuon = PFMuonAlgo::isMuon(elt);
02624 bool thisIsAnIsolatedMuon = PFMuonAlgo::isIsolatedMuon(elt);
02625 bool thisIsAGlobalTightMuon = PFMuonAlgo::isGlobalTightMuon(elt);
02626 bool thisIsATrackerTightMuon = PFMuonAlgo::isTrackerTightMuon(elt);
02627
02628
02629 bool isFromDisp = isFromSecInt(elt, "secondary");
02630 bool isToDisp = isFromSecInt(elt, "primary");
02631
02632 bool globalFitUsed = false;
02633
02634 if ( thisIsAMuon ) {
02635
02636
02637 px = muonRef->px();
02638 py = muonRef->py();
02639 pz = muonRef->pz();
02640 energy = sqrt(muonRef->p()*muonRef->p() + 0.1057*0.1057);
02641
02642 reco::TrackBase::TrackQuality trackQualityHighPurity = TrackBase::qualityByName("highPurity");
02643 if(debug_)if(!trackRef->quality(trackQualityHighPurity))cout<<" Low Purity Track "<<endl;
02644
02645 if(muonRef->isGlobalMuon() && !usePFMuonMomAssign_){
02646
02647 reco::TrackRef combinedMu = muonRef->combinedMuon();
02648
02649
02650 bool useGlobalFit = false;
02651
02652 if(thisIsAnIsolatedMuon && (!muonRef->isTrackerMuon() || (muonRef->pt() > combinedMu->pt() && track.ptError() > 5.0*combinedMu->ptError()))) useGlobalFit = true;
02653 else if(!trackRef->quality(trackQualityHighPurity)) useGlobalFit = true;
02654 else if(muonRef->pt() > combinedMu->pt() &&
02655 (track.hitPattern().numberOfValidTrackerHits() < 8 || track.hitPattern().numberOfValidPixelHits() == 0 ) &&
02656 track.ptError() > 5.0*combinedMu->ptError()) useGlobalFit = true;
02657
02658 if(useGlobalFit){
02659 px = combinedMu->px();
02660 py = combinedMu->py();
02661 pz = combinedMu->pz();
02662 energy = sqrt(combinedMu->p()*combinedMu->p() + 0.1057*0.1057);
02663 globalFitUsed = true;
02664 }
02665 }
02666 else if(usePFMuonMomAssign_){
02667
02668 if(thisIsAGlobalTightMuon)
02669 {
02670
02671 if(muonRef->isTrackerMuon()){
02672 if(sqrt(px*px+py*py) > 10){
02673 reco::TrackRef combinedMu = muonRef->combinedMuon();
02674 px = combinedMu->px();
02675 py = combinedMu->py();
02676 pz = combinedMu->pz();
02677 energy = sqrt(combinedMu->p()*combinedMu->p() + 0.1057*0.1057);
02678 globalFitUsed = true;
02679 }
02680 }
02681 else{
02682 reco::TrackRef combinedMu =
02683 muonRef->combinedMuon()->normalizedChi2() < muonRef->standAloneMuon()->normalizedChi2() ?
02684 muonRef->combinedMuon() :
02685 muonRef->standAloneMuon() ;
02686 px = combinedMu->px();
02687 py = combinedMu->py();
02688 pz = combinedMu->pz();
02689 energy = sqrt(combinedMu->p()*combinedMu->p() + 0.1057*0.1057);
02690 globalFitUsed = true;
02691 }
02692 }
02693 }
02694 }
02695 else if (isFromDisp) {
02696 double Dpt = trackRef->ptError();
02697 double dptRel = Dpt/trackRef->pt()*100;
02698
02699
02700
02701 if (dptRel < dptRel_DispVtx_){
02702
02703 if (debug_)
02704 cout << "Not refitted px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy << endl;
02705
02706 reco::PFDisplacedVertexRef vRef = eltTrack->displacedVertexRef(reco::PFBlockElement::T_FROM_DISP)->displacedVertexRef();
02707 reco::Track trackRefit = vRef->refittedTrack(trackRef);
02708 px = trackRefit.px();
02709 py = trackRefit.py();
02710 pz = trackRefit.pz();
02711 energy = sqrt(trackRefit.p()*trackRefit.p() + 0.13957*0.13957);
02712 if (debug_)
02713 cout << "Refitted px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy << endl;
02714
02715 }
02716 }
02717
02718 if ((isFromDisp || isToDisp) && debug_) cout << "Final px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy << endl;
02719
02720
02721 math::XYZTLorentzVector momentum(px,py,pz,energy);
02722 reco::PFCandidate::ParticleType particleType
02723 = reco::PFCandidate::h;
02724
02725
02726 pfCandidates_->push_back( PFCandidate( charge,
02727 momentum,
02728 particleType ) );
02729
02730
02731 if( isFromDisp ) {
02732 pfCandidates_->back().setFlag( reco::PFCandidate::T_FROM_DISP, true);
02733 pfCandidates_->back().setDisplacedVertexRef( eltTrack->displacedVertexRef(reco::PFBlockElement::T_FROM_DISP)->displacedVertexRef(), reco::PFCandidate::T_FROM_DISP);
02734 }
02735
02736
02737 if( isToDisp && !thisIsAMuon ) {
02738 pfCandidates_->back().setFlag( reco::PFCandidate::T_TO_DISP, true);
02739 pfCandidates_->back().setDisplacedVertexRef( eltTrack->displacedVertexRef(reco::PFBlockElement::T_TO_DISP)->displacedVertexRef(), reco::PFCandidate::T_TO_DISP);
02740 }
02741
02742
02743
02744 if ( thisIsAMuon && globalFitUsed )
02745 pfCandidates_->back().setVertexSource( PFCandidate::kComMuonVertex );
02746 else
02747 pfCandidates_->back().setVertexSource( PFCandidate::kTrkVertex );
02748
02749 pfCandidates_->back().setTrackRef( trackRef );
02750 pfCandidates_->back().setPositionAtECALEntrance( eltTrack->positionAtECALEntrance());
02751
02752
02753
02754
02755 if (muonRef.isNonnull()) {
02756 pfCandidates_->back().setMuonRef( muonRef );
02757
02758 if ( thisIsAMuon) {
02759 particleType = reco::PFCandidate::mu;
02760 pfCandidates_->back().setParticleType( particleType );
02761 if (debug_) {
02762 if(thisIsAGlobalTightMuon) cout << "PFAlgo: particle type set to muon (global, tight), pT = " <<muonRef->pt()<< endl;
02763 else if(thisIsATrackerTightMuon) cout << "PFAlgo: particle type set to muon (tracker, tight), pT = " <<muonRef->pt()<< endl;
02764 else if(thisIsAnIsolatedMuon) cout << "PFAlgo: particle type set to muon (isolated), pT = " <<muonRef->pt()<< endl;
02765 else cout<<" problem with muon assignment "<<endl;
02766 PFMuonAlgo::printMuonProperties( muonRef );
02767 }
02768 }
02769 }
02770
02771
02772
02773 if(debug_)
02774 cout<<"** candidate: "<<pfCandidates_->back()<<endl;
02775
02776
02777 return pfCandidates_->size()-1;
02778 }
02779
02780
02781
02782
02783
02784 unsigned
02785 PFAlgo::reconstructCluster(const reco::PFCluster& cluster,
02786 double particleEnergy,
02787 bool useDirection,
02788 double particleX,
02789 double particleY,
02790 double particleZ) {
02791
02792 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::X;
02793
02794
02795
02796
02797
02798 double factor = 1.;
02799 if ( useDirection ) {
02800 switch( cluster.layer() ) {
02801 case PFLayer::ECAL_BARREL:
02802 case PFLayer::HCAL_BARREL1:
02803 factor = std::sqrt(cluster.position().Perp2()/(particleX*particleX+particleY*particleY));
02804 break;
02805 case PFLayer::ECAL_ENDCAP:
02806 case PFLayer::HCAL_ENDCAP:
02807 case PFLayer::HF_HAD:
02808 case PFLayer::HF_EM:
02809 factor = cluster.position().Z()/particleZ;
02810 break;
02811 default:
02812 assert(0);
02813 }
02814 }
02815
02816 math::XYZPoint vertexPos;
02817 if(useVertices_)
02818 vertexPos = math::XYZPoint(primaryVertex_.x(),primaryVertex_.y(),primaryVertex_.z());
02819 else
02820 vertexPos = math::XYZPoint(0.0,0.0,0.0);
02821
02822
02823 math::XYZVector clusterPos( cluster.position().X()-vertexPos.X(),
02824 cluster.position().Y()-vertexPos.Y(),
02825 cluster.position().Z()-vertexPos.Z());
02826 math::XYZVector particleDirection ( particleX*factor-vertexPos.X(),
02827 particleY*factor-vertexPos.Y(),
02828 particleZ*factor-vertexPos.Z() );
02829
02830
02831
02832
02833 clusterPos = useDirection ? particleDirection.Unit() : clusterPos.Unit();
02834 clusterPos *= particleEnergy;
02835
02836
02837
02838
02839 double mass = 0;
02840 ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> >
02841 momentum( clusterPos.X(), clusterPos.Y(), clusterPos.Z(), mass);
02842
02843 math::XYZTLorentzVector tmp;
02844
02845 tmp = momentum;
02846
02847
02848 int charge = 0;
02849
02850
02851 switch( cluster.layer() ) {
02852 case PFLayer::ECAL_BARREL:
02853 case PFLayer::ECAL_ENDCAP:
02854 particleType = PFCandidate::gamma;
02855 break;
02856 case PFLayer::HCAL_BARREL1:
02857 case PFLayer::HCAL_ENDCAP:
02858 particleType = PFCandidate::h0;
02859 break;
02860 case PFLayer::HF_HAD:
02861 particleType = PFCandidate::h_HF;
02862 break;
02863 case PFLayer::HF_EM:
02864 particleType = PFCandidate::egamma_HF;
02865 break;
02866 default:
02867 assert(0);
02868 }
02869
02870
02871 pfCandidates_->push_back( PFCandidate( charge,
02872 tmp,
02873 particleType ) );
02874
02875
02876
02877 pfCandidates_->back().
02878 setPositionAtECALEntrance(math::XYZPointF(cluster.position().X(),
02879 cluster.position().Y(),
02880 cluster.position().Z()));
02881
02882
02883 pfCandidates_->back().setVertex(vertexPos);
02884
02885 if(debug_)
02886 cout<<"** candidate: "<<pfCandidates_->back()<<endl;
02887
02888
02889 return pfCandidates_->size()-1;
02890
02891 }
02892
02893
02894
02895
02896 double
02897 PFAlgo::neutralHadronEnergyResolution(double clusterEnergyHCAL, double eta) const {
02898
02899
02900 if ( clusterEnergyHCAL < 1. ) clusterEnergyHCAL = 1.;
02901
02902 double resol = fabs(eta) < 1.48 ?
02903 sqrt (1.02*1.02/clusterEnergyHCAL + 0.065*0.065)
02904 :
02905 sqrt (1.20*1.20/clusterEnergyHCAL + 0.028*0.028);
02906
02907 return resol;
02908
02909 }
02910
02911 double
02912 PFAlgo::nSigmaHCAL(double clusterEnergyHCAL, double eta) const {
02913 double nS = fabs(eta) < 1.48 ?
02914 nSigmaHCAL_ * (1. + exp(-clusterEnergyHCAL/100.))
02915 :
02916 nSigmaHCAL_ * (1. + exp(-clusterEnergyHCAL/100.));
02917
02918 return nS;
02919 }
02920
02921
02922 ostream& operator<<(ostream& out, const PFAlgo& algo) {
02923 if(!out ) return out;
02924
02925 out<<"====== Particle Flow Algorithm ======= ";
02926 out<<endl;
02927 out<<"nSigmaECAL_ "<<algo.nSigmaECAL_<<endl;
02928 out<<"nSigmaHCAL_ "<<algo.nSigmaHCAL_<<endl;
02929 out<<endl;
02930 out<<(*(algo.calibration_))<<endl;
02931 out<<endl;
02932 out<<"reconstructed particles: "<<endl;
02933
02934 const std::auto_ptr< reco::PFCandidateCollection >&
02935 candidates = algo.pfCandidates();
02936
02937 if(!candidates.get() ) {
02938 out<<"candidates already transfered"<<endl;
02939 return out;
02940 }
02941 for(PFCandidateConstIterator ic=algo.pfCandidates_->begin();
02942 ic != algo.pfCandidates_->end(); ic++) {
02943 out<<(*ic)<<endl;
02944 }
02945
02946 return out;
02947 }
02948
02949
02950 bool PFAlgo::isSatelliteCluster( const reco::PFRecTrack& track,
02951 const PFCluster& cluster ){
02952
02953
02954
02955
02956
02957
02958
02959 bool link = false;
02960
02961 if( cluster.layer() == PFLayer::HCAL_BARREL1 ||
02962 cluster.layer() == PFLayer::HCAL_ENDCAP ) {
02963
02964 const reco::PFTrajectoryPoint& atHCAL
02965 = track.extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance );
02966
02967 if( atHCAL.isValid() ){
02968 double tracketa = atHCAL.position().Eta();
02969 double trackphi = atHCAL.position().Phi();
02970 double hcaleta = cluster.positionREP().Eta();
02971 double hcalphi = cluster.positionREP().Phi();
02972
02973
02974 double deta = hcaleta - tracketa;
02975 double dphi = acos(cos(hcalphi - trackphi));
02976 double dr = sqrt(deta*deta + dphi*dphi);
02977
02978 if( debug_ ){
02979 cout<<"\t\t\tSatellite Test "
02980 <<tracketa<<" "<<trackphi<<" "
02981 <<hcaleta<<" "<<hcalphi<<" dr="<<dr
02982 <<endl;
02983 }
02984
02985
02986
02987
02988 if( dr < 0.17 ) link = true;
02989 }
02990
02991 }
02992
02993 return link;
02994 }
02995
02996 reco::PFBlockRef
02997 PFAlgo::createBlockRef( const reco::PFBlockCollection& blocks,
02998 unsigned bi ) {
02999
03000 if( blockHandle_.isValid() ) {
03001 return reco::PFBlockRef( blockHandle_, bi );
03002 }
03003 else {
03004 return reco::PFBlockRef( &blocks, bi );
03005 }
03006 }
03007
03008 void
03009 PFAlgo::associatePSClusters(unsigned iEcal,
03010 reco::PFBlockElement::Type psElementType,
03011 const reco::PFBlock& block,
03012 const edm::OwnVector< reco::PFBlockElement >& elements,
03013 const reco::PFBlock::LinkData& linkData,
03014 std::vector<bool>& active,
03015 std::vector<double>& psEne) {
03016
03017
03018
03019
03020
03021
03022
03023 std::multimap<double, unsigned> sortedPS;
03024 typedef std::multimap<double, unsigned>::iterator IE;
03025 block.associatedElements( iEcal, linkData,
03026 sortedPS, psElementType,
03027 reco::PFBlock::LINKTEST_ALL );
03028
03029
03030 double totalPS = 0.;
03031 for ( IE ips=sortedPS.begin(); ips!=sortedPS.end(); ++ips ) {
03032
03033
03034 unsigned iPS = ips->second;
03035
03036
03037
03038 if (!active[iPS]) continue;
03039
03040
03041 std::multimap<double, unsigned> sortedECAL;
03042 block.associatedElements( iPS, linkData,
03043 sortedECAL,
03044 reco::PFBlockElement::ECAL,
03045 reco::PFBlock::LINKTEST_ALL );
03046 unsigned jEcal = sortedECAL.begin()->second;
03047 if ( jEcal != iEcal ) continue;
03048
03049
03050 PFBlockElement::Type pstype = elements[ iPS ].type();
03051 assert( pstype == psElementType );
03052 PFClusterRef psclusterref = elements[iPS].clusterRef();
03053 assert(!psclusterref.isNull() );
03054 totalPS += psclusterref->energy();
03055 psEne[0] += psclusterref->energy();
03056 active[iPS] = false;
03057 }
03058
03059
03060 }
03061
03062
03063 bool
03064 PFAlgo::isFromSecInt(const reco::PFBlockElement& eTrack, string order) const {
03065
03066 reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
03067 reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
03068 reco::PFBlockElement::TrackType T_FROM_GAMMACONV = reco::PFBlockElement::T_FROM_GAMMACONV;
03069 reco::PFBlockElement::TrackType T_FROM_V0 = reco::PFBlockElement::T_FROM_V0;
03070
03071 bool bPrimary = (order.find("primary") != string::npos);
03072 bool bSecondary = (order.find("secondary") != string::npos);
03073 bool bAll = (order.find("all") != string::npos);
03074
03075 bool isToDisp = usePFNuclearInteractions_ && eTrack.trackType(T_TO_DISP);
03076 bool isFromDisp = usePFNuclearInteractions_ && eTrack.trackType(T_FROM_DISP);
03077
03078 if (bPrimary && isToDisp) return true;
03079 if (bSecondary && isFromDisp ) return true;
03080 if (bAll && ( isToDisp || isFromDisp ) ) return true;
03081
03082 bool isFromConv = usePFConversions_ && eTrack.trackType(T_FROM_GAMMACONV);
03083
03084 if ((bAll || bSecondary)&& isFromConv) return true;
03085
03086 bool isFromDecay = (bAll || bSecondary) && usePFDecays_ && eTrack.trackType(T_FROM_V0);
03087
03088 return isFromDecay;
03089
03090
03091 }
03092
03093
03094 void
03095 PFAlgo::setEGElectronCollection(const reco::GsfElectronCollection & egelectrons) {
03096 if(useEGElectrons_ && pfele_) pfele_->setEGElectronCollection(egelectrons);
03097 }
03098
03099 void
03100 PFAlgo::postCleaning() {
03101
03102
03103 if ( !postHFCleaning_ ) return;
03104
03105
03106 double metX = 0.;
03107 double metY = 0.;
03108 double sumet = 0;
03109 std::vector<unsigned int> pfCandidatesToBeRemoved;
03110 for(unsigned i=0; i<pfCandidates_->size(); i++) {
03111 const PFCandidate& pfc = (*pfCandidates_)[i];
03112 metX += pfc.px();
03113 metY += pfc.py();
03114 sumet += pfc.pt();
03115 }
03116 double met2 = metX*metX+metY*metY;
03117
03118 double significance = std::sqrt(met2/sumet);
03119 double significanceCor = significance;
03120 if ( significance > minSignificance_ ) {
03121
03122 double metXCor = metX;
03123 double metYCor = metY;
03124 double sumetCor = sumet;
03125 double met2Cor = met2;
03126 double deltaPhi = 3.14159;
03127 double deltaPhiPt = 100.;
03128 double deltaPhiCor = 3.14159;
03129 double deltaPhiPtCor = 100.;
03130 bool next = true;
03131 unsigned iCor = 1E9;
03132
03133
03134 while ( next ) {
03135
03136 double metReduc = -1.;
03137 double setReduc = -1.;
03138
03139 for(unsigned i=0; i<pfCandidates_->size(); ++i) {
03140 const PFCandidate& pfc = (*pfCandidates_)[i];
03141
03142
03143 if ( pfc.particleId() != reco::PFCandidate::h_HF &&
03144 pfc.particleId() != reco::PFCandidate::egamma_HF ) continue;
03145
03146
03147 if ( pfc.pt() < minHFCleaningPt_ ) continue;
03148
03149
03150 bool skip = false;
03151 for(unsigned j=0; j<pfCandidatesToBeRemoved.size(); ++j) {
03152 if ( i == pfCandidatesToBeRemoved[j] ) skip = true;
03153 if ( skip ) break;
03154 }
03155 if ( skip ) continue;
03156
03157
03158 deltaPhi = std::acos((metX*pfc.px()+metY*pfc.py())/(pfc.pt()*std::sqrt(met2)));
03159 deltaPhiPt = deltaPhi*pfc.pt();
03160 if ( deltaPhiPt > maxDeltaPhiPt_ ) continue;
03161
03162
03163 double metXInt = metX - pfc.px();
03164 double metYInt = metY - pfc.py();
03165 double sumetInt = sumet - pfc.pt();
03166 double met2Int = metXInt*metXInt+metYInt*metYInt;
03167 if ( met2Int < met2Cor ) {
03168 metXCor = metXInt;
03169 metYCor = metYInt;
03170 metReduc = (met2-met2Int)/met2Int;
03171 setReduc = (std::sqrt(met2Int)-std::sqrt(met2))/(sumetInt-sumet);
03172 met2Cor = met2Int;
03173 sumetCor = sumetInt;
03174 significanceCor = std::sqrt(met2Cor/sumetCor);
03175 iCor = i;
03176 deltaPhiCor = deltaPhi;
03177 deltaPhiPtCor = deltaPhiPt;
03178 }
03179 }
03180
03181
03182 if ( metReduc > minDeltaMet_ ) {
03183 pfCandidatesToBeRemoved.push_back(iCor);
03184 metX = metXCor;
03185 metY = metYCor;
03186 sumet = sumetCor;
03187 met2 = met2Cor;
03188 } else {
03189
03190 next = false;
03191 }
03192 }
03193
03194
03195 if ( significance - significanceCor > minSignificanceReduction_ &&
03196 significanceCor < maxSignificance_ ) {
03197 std::cout << "Significance reduction = " << significance << " -> "
03198 << significanceCor << " = " << significanceCor - significance
03199 << std::endl;
03200 for(unsigned j=0; j<pfCandidatesToBeRemoved.size(); ++j) {
03201 std::cout << "Removed : " << (*pfCandidates_)[pfCandidatesToBeRemoved[j]] << std::endl;
03202 pfCleanedCandidates_->push_back( (*pfCandidates_)[ pfCandidatesToBeRemoved[j] ] );
03203 (*pfCandidates_)[pfCandidatesToBeRemoved[j]].rescaleMomentum(1E-6);
03204
03205
03206 }
03207 }
03208 }
03209
03210 }
03211
03212 void
03213 PFAlgo::checkCleaning( const reco::PFRecHitCollection& cleanedHits ) {
03214
03215
03216
03217 if ( !cleanedHits.size() ) return;
03218
03219
03220 double metX = 0.;
03221 double metY = 0.;
03222 double sumet = 0;
03223 std::vector<unsigned int> hitsToBeAdded;
03224 for(unsigned i=0; i<pfCandidates_->size(); i++) {
03225 const PFCandidate& pfc = (*pfCandidates_)[i];
03226 metX += pfc.px();
03227 metY += pfc.py();
03228 sumet += pfc.pt();
03229 }
03230 double met2 = metX*metX+metY*metY;
03231 double met2_Original = met2;
03232
03233
03234
03235 double metXCor = metX;
03236 double metYCor = metY;
03237 double sumetCor = sumet;
03238 double met2Cor = met2;
03239 bool next = true;
03240 unsigned iCor = 1E9;
03241
03242
03243 while ( next ) {
03244
03245 double metReduc = -1.;
03246 double setReduc = -1.;
03247
03248 for(unsigned i=0; i<cleanedHits.size(); ++i) {
03249 const PFRecHit& hit = cleanedHits[i];
03250 double length = std::sqrt(hit.position().Mag2());
03251 double px = hit.energy() * hit.position().x()/length;
03252 double py = hit.energy() * hit.position().y()/length;
03253 double pt = std::sqrt(px*px + py*py);
03254
03255
03256 bool skip = false;
03257 for(unsigned j=0; j<hitsToBeAdded.size(); ++j) {
03258 if ( i == hitsToBeAdded[j] ) skip = true;
03259 if ( skip ) break;
03260 }
03261 if ( skip ) continue;
03262
03263
03264 double metXInt = metX + px;
03265 double metYInt = metY + py;
03266 double sumetInt = sumet + pt;
03267 double met2Int = metXInt*metXInt+metYInt*metYInt;
03268
03269
03270 if ( met2Int < met2Cor ) {
03271 metXCor = metXInt;
03272 metYCor = metYInt;
03273 metReduc = (met2-met2Int)/met2Int;
03274 setReduc = (std::sqrt(met2Int)-std::sqrt(met2))/(sumetInt-sumet);
03275 met2Cor = met2Int;
03276 sumetCor = sumetInt;
03277
03278 iCor = i;
03279 }
03280 }
03281
03282
03283
03284 if ( metReduc > minDeltaMet_ ) {
03285 hitsToBeAdded.push_back(iCor);
03286 metX = metXCor;
03287 metY = metYCor;
03288 sumet = sumetCor;
03289 met2 = met2Cor;
03290 } else {
03291
03292 next = false;
03293 }
03294 }
03295
03296
03297 if ( std::sqrt(met2_Original) - std::sqrt(met2) > 5. ) {
03298 if ( debug_ ) {
03299 std::cout << hitsToBeAdded.size() << " hits were re-added " << std::endl;
03300 std::cout << "MET reduction = " << std::sqrt(met2_Original) << " -> "
03301 << std::sqrt(met2Cor) << " = " << std::sqrt(met2Cor) - std::sqrt(met2_Original)
03302 << std::endl;
03303 std::cout << "Added after cleaning check : " << std::endl;
03304 }
03305 for(unsigned j=0; j<hitsToBeAdded.size(); ++j) {
03306 const PFRecHit& hit = cleanedHits[hitsToBeAdded[j]];
03307 PFCluster cluster(hit.layer(), hit.energy(),
03308 hit.position().x(), hit.position().y(), hit.position().z() );
03309 reconstructCluster(cluster,hit.energy());
03310 if ( debug_ ) {
03311 std::cout << pfCandidates_->back() << ". time = " << hit.time() << std::endl;
03312 }
03313 }
03314 }
03315
03316 }
03317
03318
03319 void
03320 PFAlgo::postMuonCleaning( const edm::Handle<reco::MuonCollection>& muonh,
03321 const reco::VertexCollection& primaryVertices) {
03322
03323
03324 bool printout = false;
03325
03326
03327
03328
03329 double sumetPU = 0.;
03330 for (unsigned short i=1 ;i<primaryVertices.size();++i ) {
03331 if ( !primaryVertices[i].isValid() || primaryVertices[i].isFake() ) continue;
03332 primaryVertices[i];
03333 for ( reco::Vertex::trackRef_iterator itr = primaryVertices[i].tracks_begin();
03334 itr < primaryVertices[i].tracks_end(); ++itr ) {
03335 sumetPU += (*itr)->pt();
03336 }
03337 }
03338 sumetPU /= 0.65;
03339
03340
03341
03342 double metX = 0.;
03343 double metY = 0.;
03344 double sumet = 0;
03345 double metXCosmics = 0.;
03346 double metYCosmics = 0.;
03347 double sumetCosmics = 0.;
03348
03349 std::map<double,unsigned int> pfMuons;
03350 std::map<double,unsigned int> pfCosmics;
03351 typedef std::multimap<double, unsigned int>::iterator IE;
03352
03353 for(unsigned i=0; i<pfCandidates_->size(); i++) {
03354 const PFCandidate& pfc = (*pfCandidates_)[i];
03355 double origin =
03356 (pfc.vx()-primaryVertex_.x())*(pfc.vx()-primaryVertex_.x()) +
03357 (pfc.vy()-primaryVertex_.y())*(pfc.vy()-primaryVertex_.y());
03358
03359
03360 metX += pfc.px();
03361 metY += pfc.py();
03362 sumet += pfc.pt();
03363
03364
03365
03366 if ( pfc.particleId() == reco::PFCandidate::mu ) {
03367 pfMuons.insert(std::pair<double,unsigned int>(-pfc.pt(),i));
03368 if ( origin > 1.0 ) {
03369 pfCosmics.insert(std::pair<double,unsigned int>(-pfc.pt(),i));
03370 metXCosmics += pfc.px();
03371 metYCosmics += pfc.py();
03372 sumetCosmics += pfc.pt();
03373 }
03374 }
03375
03376 }
03377
03378 double met2 = metX*metX+metY*metY;
03379 double met2Cosmics = (metX-metXCosmics)*(metX-metXCosmics)+(metY-metYCosmics)*(metY-metYCosmics);
03380
03381
03382 if ( sumetCosmics > (sumet-sumetPU)/10. && met2Cosmics < met2 ) {
03383 if ( printout )
03384 std::cout << "MEX,MEY,MET Before (Cosmics)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03385 for ( IE imu = pfCosmics.begin(); imu != pfCosmics.end(); ++imu ) {
03386 const PFCandidate& pfc = (*pfCandidates_)[imu->second];
03387
03388 pfCosmicsMuonCleanedCandidates_->push_back(pfc);
03389 if ( printout )
03390 std::cout << "Muon cleaned (cosmic). pt/d0 = " << pfc.pt() << " "
03391 << std::sqrt(pfc.vx()*pfc.vx() + pfc.vy()*pfc.vy()) << std::endl;
03392 metX -= pfc.px();
03393 metY -= pfc.py();
03394 sumet -= pfc.pt();
03395 (*pfCandidates_)[imu->second].rescaleMomentum(1E-6);
03396 }
03397 met2 = metX*metX+metY*metY;
03398 if ( printout )
03399 std::cout << "MEX,MEY,MET After (Cosmics)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03400 }
03401
03402
03403
03404
03405
03406 for ( IE imu = pfMuons.begin(); imu != pfMuons.end(); ++imu ) {
03407 const PFCandidate& pfc = (*pfCandidates_)[imu->second];
03408
03409
03410 if ( pfc.pt() < 20. ) continue;
03411
03412 double metXNO = metX - pfc.px();
03413 double metYNO = metY - pfc.py();
03414 double met2NO = metXNO*metXNO + metYNO*metYNO;
03415 double sumetNO = sumet - pfc.pt();
03416 double factor = std::max(2.,2000./sumetNO);
03417
03418 reco::MuonRef muonRef = pfc.muonRef();
03419 reco::TrackRef combinedMu = muonRef->combinedMuon();
03420 reco::TrackRef trackerMu = muonRef->track();
03421 reco::TrackRef standAloneMu = muonRef->standAloneMuon();
03422
03423 if ( !combinedMu || !trackerMu ) {
03424 if ( sumetNO-sumetPU > 500. && met2NO < met2/4.) {
03425 pfFakeMuonCleanedCandidates_->push_back(pfc);
03426 if ( printout ) {
03427 std::cout << "Muon cleaned (NO-bad) " << sumetNO << std::endl;
03428 std::cout << "sumet NO " << sumetNO << std::endl;
03429 }
03430 (*pfCandidates_)[imu->second].rescaleMomentum(1E-6);
03431 if ( printout )
03432 std::cout << "MEX,MEY,MET " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03433 metX = metXNO;
03434 metY = metYNO;
03435 met2 = met2NO;
03436 sumet = sumetNO;
03437 if ( printout ) {
03438 std::cout << "Muon cleaned (NO/TK) ! " << std::endl;
03439 std::cout << "MEX,MEY,MET Now (NO/TK)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03440 }
03441 }
03442 } else {
03443
03444
03445
03446
03447
03448
03449
03450
03451 double metXTK = metXNO + trackerMu->px();
03452 double metYTK = metYNO + trackerMu->py();
03453 double met2TK = metXTK*metXTK + metYTK*metYTK;
03454
03455 double metXGL = metXNO + combinedMu->px();
03456 double metYGL = metYNO + combinedMu->py();
03457 double met2GL = metXGL*metXGL + metYGL*metYGL;
03458
03459
03460
03461
03462
03463
03464
03465
03466
03467
03468
03469
03470
03471
03472
03473 bool fixed = false;
03474 if ( ( sumetNO-sumetPU > 250. && met2TK < met2/4. && met2TK < met2GL ) ||
03475 ( met2TK < met2/2. && trackerMu->pt() < combinedMu->pt()/4. && met2TK < met2GL ) ) {
03476 pfCleanedTrackerAndGlobalMuonCandidates_->push_back(pfc);
03477 math::XYZTLorentzVectorD p4(trackerMu->px(),trackerMu->py(),trackerMu->pz(),
03478 std::sqrt(trackerMu->p()*trackerMu->p()+0.1057*0.1057));
03479 if ( printout )
03480 std::cout << "Muon cleaned (TK) ! " << met2TK/met2 << " " << trackerMu->pt() / combinedMu->pt() << std::endl;
03481 (*pfCandidates_)[imu->second].setP4(p4);
03482 if ( printout )
03483 std::cout << "MEX,MEY,MET Before (TK)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03484 metX = metXTK;
03485 metY = metYTK;
03486 met2 = met2TK;
03487 fixed = true;
03488 if ( printout )
03489 std::cout << "MEX,MEY,MET Now (TK)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03490 }
03491
03492 else if ( ( sumetNO-sumetPU > 250. && met2GL < met2/4. && met2GL < met2TK ) ||
03493 ( met2GL < met2/2. && combinedMu->pt() < trackerMu->pt()/4.&& met2GL < met2TK ) ) {
03494 pfCleanedTrackerAndGlobalMuonCandidates_->push_back(pfc);
03495 math::XYZTLorentzVectorD p4(combinedMu->px(),combinedMu->py(),combinedMu->pz(),
03496 std::sqrt(combinedMu->p()*combinedMu->p()+0.1057*0.1057));
03497 if ( printout )
03498 std::cout << "Muon cleaned (GL) ! " << met2GL/met2 << " " << combinedMu->pt()/trackerMu->pt() << std::endl;
03499 (*pfCandidates_)[imu->second].setP4(p4);
03500 if ( printout )
03501 std::cout << "MEX,MEY,MET before (GL)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03502 metX = metXGL;
03503 metY = metYGL;
03504 met2 = met2GL;
03505 fixed = true;
03506 if ( printout )
03507 std::cout << "MEX,MEY,MET Now (GL)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03508 }
03509
03510
03511 bool fake =
03512 fabs ( pfc.eta() ) > 2.15 &&
03513 met2NO < met2/25. &&
03514 (met2GL < met2TK/2. || met2TK < met2GL/2.) &&
03515 standAloneMu->pt() < combinedMu->pt()/10. &&
03516 standAloneMu->pt() < trackerMu->pt()/10.;
03517
03518
03519 bool punchthrough1 =
03520 ( sumetNO-sumetPU > 250. && met2NO < met2/4. && (met2GL < met2TK/factor || met2TK < met2GL/factor) );
03521
03522
03523 bool punchthrough2 =
03524 pfc.p() > 100. && pfc.rawHcalEnergy() > 100. &&
03525 pfc.rawEcalEnergy()+pfc.rawHcalEnergy() > pfc.p()/3. &&
03526 !PFMuonAlgo::isIsolatedMuon(muonRef) && met2NO < met2/4.;
03527
03528 if ( punchthrough1 || punchthrough2 || fake ) {
03529
03530
03531 const PFCandidate::ElementsInBlocks& eleInBlocks = pfc.elementsInBlocks();
03532 if ( !eleInBlocks.size() ) {
03533 ostringstream err;
03534 err<<"A muon has no associated elements in block. Cannot proceed with postMuonCleaning()";
03535 edm::LogError("PFAlgo")<<err.str()<<endl;
03536 } else {
03537 PFBlockRef blockRefMuon = eleInBlocks[0].first;
03538 unsigned indexMuon = eleInBlocks[0].second;
03539 for ( unsigned iele = 1; iele < eleInBlocks.size(); ++iele ) {
03540 indexMuon = eleInBlocks[iele].second;
03541 break;
03542 }
03543
03544
03545 double iHad = 1E9;
03546 bool hadron = false;
03547 for ( unsigned i = imu->second+1; i < pfCandidates_->size(); ++i ) {
03548 const PFCandidate& pfcn = (*pfCandidates_)[i];
03549 const PFCandidate::ElementsInBlocks& ele = pfcn.elementsInBlocks();
03550 if ( !ele.size() ) {
03551 ostringstream err2;
03552 err2<<"A pfCandidate has no associated elements in block. Cannot proceed with postMuonCleaning()";
03553 edm::LogError("PFAlgo")<<err2.str()<<endl;
03554 continue;
03555 }
03556 PFBlockRef blockRefHadron = ele[0].first;
03557 unsigned indexHadron = ele[0].second;
03558
03559 if ( blockRefHadron.key() != blockRefMuon.key() ) break;
03560
03561 if ( indexHadron == indexMuon &&
03562 pfcn.particleId() == reco::PFCandidate::h0 ) {
03563 iHad = i;
03564 hadron = true;
03565 }
03566 if ( hadron ) break;
03567 }
03568
03569 if ( hadron ) {
03570 const PFCandidate& pfch = (*pfCandidates_)[iHad];
03571 pfPunchThroughMuonCleanedCandidates_->push_back(pfc);
03572 pfPunchThroughHadronCleanedCandidates_->push_back(pfch);
03573 if ( printout ) {
03574 std::cout << pfc << std::endl;
03575 std::cout << "Raw ecal/hcal : " << pfc.rawEcalEnergy() << " " << pfc.rawHcalEnergy() << std::endl;
03576 std::cout << "Hadron: " << (*pfCandidates_)[iHad] << std::endl;
03577 }
03578 double rescaleFactor = (*pfCandidates_)[iHad].p()/(*pfCandidates_)[imu->second].p();
03579 metX -= (*pfCandidates_)[imu->second].px() + (*pfCandidates_)[iHad].px();
03580 metY -= (*pfCandidates_)[imu->second].py() + (*pfCandidates_)[iHad].py();
03581 (*pfCandidates_)[imu->second].rescaleMomentum(rescaleFactor);
03582 (*pfCandidates_)[iHad].rescaleMomentum(1E-6);
03583 (*pfCandidates_)[imu->second].setParticleType(reco::PFCandidate::h);
03584 if ( printout )
03585 std::cout << "MEX,MEY,MET Before " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03586 metX += (*pfCandidates_)[imu->second].px();
03587 metY += (*pfCandidates_)[imu->second].py();
03588 met2 = metX*metX + metY*metY;
03589 if ( printout ) {
03590 std::cout << "Muon changed to charged hadron" << std::endl;
03591 std::cout << "MEX,MEY,MET Now (NO)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03592 }
03593 } else if ( punchthrough1 || fake ) {
03594 const PFCandidate& pfc = (*pfCandidates_)[imu->second];
03595 pfFakeMuonCleanedCandidates_->push_back(pfc);
03596 if ( printout )
03597 std::cout << "MEX,MEY,MET Before " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03598 metX -= (*pfCandidates_)[imu->second].px();
03599 metY -= (*pfCandidates_)[imu->second].py();
03600 met2 = metX*metX + metY*metY;
03601 (*pfCandidates_)[imu->second].rescaleMomentum(1E-6);
03602 if ( printout ) {
03603 std::cout << "Muon cleaned (NO)" << std::endl;
03604 std::cout << "MEX,MEY,MET Now (NO)" << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03605 }
03606 }
03607 }
03608 }
03609
03610
03611
03612
03613
03614
03615
03616
03617
03618
03619 }
03620
03621 }
03622
03623
03624 for ( unsigned imu = 0; imu < muonh->size(); ++imu ) {
03625 reco::MuonRef muonRef( muonh, imu );
03626
03627 reco::TrackRef combinedMu = muonRef->combinedMuon();
03628 reco::TrackRef trackerMu = muonRef->track();
03629 reco::TrackRef standAloneMu = muonRef->standAloneMuon();
03630
03631
03632 bool used = false;
03633 bool hadron = false;
03634 unsigned iHad = 1E9;
03635 for(unsigned i=0; i<pfCandidates_->size(); i++) {
03636 const PFCandidate& pfc = (*pfCandidates_)[i];
03637 if ( !pfc.trackRef().isNonnull() ) continue;
03638
03639 if ( pfc.trackRef().isNonnull() && pfc.trackRef() == trackerMu ) {
03640 hadron = true;
03641 iHad = i;
03642 }
03643
03644
03645
03646 if ( !pfc.muonRef().isNonnull() ) continue;
03647
03648
03649 if ( pfc.muonRef()->track() == trackerMu || pfc.muonRef()->combinedMuon() == combinedMu ) {
03650 if ( printout ) {
03651 std::cout << "This muon is already used ..." << std::endl;
03652 std::cout << pfc << std::endl;
03653 std::cout << muonRef->p() << " " << muonRef->pt() << " " << muonRef->eta() << " " << muonRef->phi() << std::endl;
03654 }
03655 used = true;
03656 }
03657 else {
03658
03659
03660 if ( pfc.muonRef()->isStandAloneMuon() && muonRef->isStandAloneMuon() ) {
03661 double dEta = pfc.muonRef()->standAloneMuon()->eta() - standAloneMu->eta();
03662 double dPhi = pfc.muonRef()->standAloneMuon()->phi() - standAloneMu->phi();
03663 double dR = sqrt(dEta*dEta + dPhi*dPhi);
03664 if ( printout ) {
03665 std::cout << "StandAlone to be added ? " << std::endl;
03666 std::cout << " eta = " << pfc.muonRef()->standAloneMuon()->eta() << " " << standAloneMu->eta() << std::endl;
03667 std::cout << " phi = " << pfc.muonRef()->standAloneMuon()->phi() << " " << standAloneMu->phi() << std::endl;
03668 std::cout << " pt = " << pfc.muonRef()->standAloneMuon()->pt() << " " << standAloneMu->pt() << std::endl;
03669 std::cout << " Delta R = " << dR << std::endl;
03670 }
03671 if ( dR < 0.005 ) {
03672 used = true;
03673 if ( printout ) std::cout << "Not removed !" << std::endl;
03674 }
03675 }
03676 }
03677 if ( used ) break;
03678
03679 }
03680
03681 if ( used ) continue;
03682
03683 double ptGL = muonRef->isGlobalMuon() ? combinedMu->pt() : 0.;
03684 double pxGL = muonRef->isGlobalMuon() ? combinedMu->px() : 0.;
03685 double pyGL = muonRef->isGlobalMuon() ? combinedMu->py() : 0.;
03686 double pzGL = muonRef->isGlobalMuon() ? combinedMu->pz() : 0.;
03687
03688 double ptTK = muonRef->isTrackerMuon() ? trackerMu->pt() : 0.;
03689 double pxTK = muonRef->isTrackerMuon() ? trackerMu->px() : 0.;
03690 double pyTK = muonRef->isTrackerMuon() ? trackerMu->py() : 0.;
03691 double pzTK = muonRef->isTrackerMuon() ? trackerMu->pz() : 0.;
03692
03693 double ptST = muonRef->isStandAloneMuon() ? standAloneMu->pt() : 0.;
03694 double pxST = muonRef->isStandAloneMuon() ? standAloneMu->px() : 0.;
03695 double pyST = muonRef->isStandAloneMuon() ? standAloneMu->py() : 0.;
03696 double pzST = muonRef->isStandAloneMuon() ? standAloneMu->pz() : 0.;
03697
03698
03699
03700 double metXTK = metX + pxTK;
03701 double metYTK = metY + pyTK;
03702 double met2TK = metXTK*metXTK + metYTK*metYTK;
03703
03704 double metXGL = metX + pxGL;
03705 double metYGL = metY + pyGL;
03706 double met2GL = metXGL*metXGL + metYGL*metYGL;
03707
03708 double metXST = metX + pxST;
03709 double metYST = metY + pyST;
03710 double met2ST = metXST*metXST + metYST*metYST;
03711
03712
03713
03714
03715 if ( ptTK > 20. && met2TK < met2/4. && met2TK < met2GL && met2TK < met2ST ) {
03716 double energy = std::sqrt(pxTK*pxTK+pyTK*pyTK+pzTK*pzTK+0.1057*0.1057);
03717 int charge = trackerMu->charge()>0 ? 1 : -1;
03718 math::XYZTLorentzVector momentum(pxTK,pyTK,pzTK,energy);
03719 reco::PFCandidate::ParticleType particleType =
03720 ptGL > 20. ? reco::PFCandidate::mu : reco::PFCandidate::h;
03721 double radius = std::sqrt(
03722 (trackerMu->vertex().x()-primaryVertex_.x())*(trackerMu->vertex().x()-primaryVertex_.x())+
03723 (trackerMu->vertex().y()-primaryVertex_.y())*(trackerMu->vertex().y()-primaryVertex_.y()));
03724
03725
03726 if ( !hadron && radius < 1.0 ) {
03727 pfCandidates_->push_back( PFCandidate( charge,
03728 momentum,
03729 particleType ) );
03730
03731 pfCandidates_->back().setVertexSource( PFCandidate::kTrkMuonVertex );
03732 pfCandidates_->back().setTrackRef( trackerMu );
03733 pfCandidates_->back().setMuonRef( muonRef );
03734
03735 if ( printout ) {
03736 std::cout << "MEX,MEY,MET before " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03737 std::cout << "Muon TK added " << std::endl;
03738 std::cout << "pT TK/GL/ST : " << ptTK << " " << ptGL << " " << ptST << std::endl;
03739 }
03740 metX += pfCandidates_->back().px();
03741 metY += pfCandidates_->back().py();
03742 met2 = metX*metX + metY*metY;
03743 if ( printout ) {
03744 std::cout << pfCandidates_->back() << std::endl;
03745 std::cout << "MEX,MEY,MET now " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03746 }
03747
03748 const PFCandidate& pfc = pfCandidates_->back();
03749 pfAddedMuonCandidates_->push_back(pfc);
03750
03751 }
03752
03753 } else if ( ptGL > 20. && met2GL < met2/4. && met2GL < met2TK && met2GL < met2ST ) {
03754
03755 double energy = std::sqrt(pxGL*pxGL+pyGL*pyGL+pzGL*pzGL+0.1057*0.1057);
03756 int charge = combinedMu->charge()>0 ? 1 : -1;
03757 math::XYZTLorentzVector momentum(pxGL,pyGL,pzGL,energy);
03758 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::mu;
03759 double radius = std::sqrt(
03760 (combinedMu->vertex().x()-primaryVertex_.x())*(combinedMu->vertex().x()-primaryVertex_.x())+
03761 (combinedMu->vertex().y()-primaryVertex_.y())*(combinedMu->vertex().y()-primaryVertex_.y()));
03762
03763
03764 if ( radius < 1.0 ) {
03765 pfCandidates_->push_back( PFCandidate( charge,
03766 momentum,
03767 particleType ) );
03768
03769 pfCandidates_->back().setVertexSource( PFCandidate::kComMuonVertex );
03770
03771 if (trackerMu.isNonnull() ) pfCandidates_->back().setTrackRef( trackerMu );
03772 pfCandidates_->back().setMuonRef( muonRef );
03773
03774 if ( printout ) {
03775 std::cout << "MEX,MEY,MET before " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03776 std::cout << "Muon GL added " << std::endl;
03777 std::cout << "pT TK/GL/ST : " << ptTK << " " << ptGL << " " << ptST << std::endl;
03778 }
03779
03780 metX += pfCandidates_->back().px();
03781 metY += pfCandidates_->back().py();
03782 met2 = metX*metX + metY*metY;
03783
03784 if ( printout ) {
03785 std::cout << pfCandidates_->back() << std::endl;
03786 std::cout << "MEX,MEY,MET now " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03787 }
03788
03789 const PFCandidate& pfc = pfCandidates_->back();
03790 pfAddedMuonCandidates_->push_back(pfc);
03791
03792 }
03793
03794 } else if ( ptST > 20. && met2ST < met2/4. && met2ST < met2TK && met2ST < met2GL ) {
03795
03796 double energy = std::sqrt(pxST*pxST+pyST*pyST+pzST*pzST+0.1057*0.1057);
03797 int charge = standAloneMu->charge()>0 ? 1 : -1;
03798 math::XYZTLorentzVector momentum(pxST,pyST,pzST,energy);
03799 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::mu;
03800 double radius = std::sqrt(
03801 (standAloneMu->vertex().x()-primaryVertex_.x())*(standAloneMu->vertex().x()-primaryVertex_.x())+
03802 (standAloneMu->vertex().y()-primaryVertex_.y())*(standAloneMu->vertex().y()-primaryVertex_.y()));
03803
03804
03805 if ( radius < 1.0 ) {
03806 pfCandidates_->push_back( PFCandidate( charge,
03807 momentum,
03808 particleType ) );
03809
03810 pfCandidates_->back().setVertexSource( PFCandidate::kSAMuonVertex);
03811 if (trackerMu.isNonnull() ) pfCandidates_->back().setTrackRef( trackerMu );
03812 pfCandidates_->back().setMuonRef( muonRef );
03813
03814 if ( printout ) {
03815 std::cout << "MEX,MEY,MET before " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03816 std::cout << "Muon ST added " << std::endl;
03817 std::cout << "pT TK/GL/ST : " << ptTK << " " << ptGL << " " << ptST << std::endl;
03818 }
03819
03820 metX += pfCandidates_->back().px();
03821 metY += pfCandidates_->back().py();
03822 met2 = metX*metX + metY*metY;
03823
03824 if ( printout ) {
03825 std::cout << pfCandidates_->back() << std::endl;
03826 std::cout << "MEX,MEY,MET now " << metX << " " << metY << " " << std::sqrt(met2) << std::endl;
03827 }
03828
03829 const PFCandidate& pfc = pfCandidates_->back();
03830 pfAddedMuonCandidates_->push_back(pfc);
03831
03832 }
03833 }
03834
03835 }
03836
03837
03838
03839
03840
03841
03842
03843
03844 }
03845
03846 void PFAlgo::setElectronExtraRef(const edm::OrphanHandle<reco::PFCandidateElectronExtraCollection >& extrah) {
03847 if(!usePFElectrons_) return;
03848
03849 unsigned size=pfCandidates_->size();
03850
03851 for(unsigned ic=0;ic<size;++ic) {
03852
03853 if((*pfCandidates_)[ic].particleId()==PFCandidate::e) {
03854
03855 PFElectronExtraEqual myExtraEqual((*pfCandidates_)[ic].gsfTrackRef());
03856 std::vector<PFCandidateElectronExtra>::const_iterator it=find_if(pfElectronExtra_.begin(),pfElectronExtra_.end(),myExtraEqual);
03857 if(it!=pfElectronExtra_.end()) {
03858
03859 reco::PFCandidateElectronExtraRef theRef(extrah,it-pfElectronExtra_.begin());
03860 (*pfCandidates_)[ic].setPFElectronExtraRef(theRef);
03861 }
03862 else {
03863 (*pfCandidates_)[ic].setPFElectronExtraRef(PFCandidateElectronExtraRef());
03864 }
03865 }
03866 else
03867 {
03868 if((*pfCandidates_)[ic].trackRef().isNonnull()) {
03869 PFElectronExtraKfEqual myExtraEqual((*pfCandidates_)[ic].trackRef());
03870 std::vector<PFCandidateElectronExtra>::const_iterator it=find_if(pfElectronExtra_.begin(),pfElectronExtra_.end(),myExtraEqual);
03871 if(it!=pfElectronExtra_.end()) {
03872 (*pfCandidates_)[ic].set_mva_e_pi(it->mvaVariable(PFCandidateElectronExtra::MVA_MVA));
03873 reco::PFCandidateElectronExtraRef theRef(extrah,it-pfElectronExtra_.begin());
03874 (*pfCandidates_)[ic].setPFElectronExtraRef(theRef);
03875 (*pfCandidates_)[ic].setGsfTrackRef(it->gsfTrackRef());
03876 }
03877 }
03878 }
03879
03880 }
03881
03882 size=pfElectronCandidates_->size();
03883 for(unsigned ic=0;ic<size;++ic) {
03884
03885 if((*pfElectronCandidates_)[ic].particleId()==PFCandidate::e) {
03886
03887 PFElectronExtraEqual myExtraEqual((*pfElectronCandidates_)[ic].gsfTrackRef());
03888 std::vector<PFCandidateElectronExtra>::const_iterator it=find_if(pfElectronExtra_.begin(),pfElectronExtra_.end(),myExtraEqual);
03889 if(it!=pfElectronExtra_.end()) {
03890 reco::PFCandidateElectronExtraRef theRef(extrah,it-pfElectronExtra_.begin());
03891 (*pfElectronCandidates_)[ic].setPFElectronExtraRef(theRef);
03892
03893 }
03894 }
03895 }
03896
03897 }