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