File indexing completed on 2024-04-06 12:01:07
0001 #include "CommonTools/ParticleFlow/interface/PFIsoDepositAlgo.h"
0002
0003 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0004
0005 #include "DataFormats/Math/interface/LorentzVector.h"
0006
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008
0009 using namespace std;
0010 using namespace edm;
0011 using namespace reco;
0012 using namespace math;
0013 using namespace pf2pat;
0014
0015 PFIsoDepositAlgo::PFIsoDepositAlgo(const edm::ParameterSet& iConfig)
0016 : verbose_(iConfig.getUntrackedParameter<bool>("verbose", false))
0017
0018 {}
0019
0020 PFIsoDepositAlgo::~PFIsoDepositAlgo() {}
0021
0022 const PFIsoDepositAlgo::IsoDeposits& PFIsoDepositAlgo::produce(const ParticleCollection& toBeIsolated,
0023 const ParticleCollection& forIsolation) {
0024 isoDeposits_.clear();
0025 isoDeposits_.reserve(toBeIsolated.size());
0026
0027 for (unsigned i = 0; i < toBeIsolated.size(); i++) {
0028 const reco::PFCandidate& toBeIso = toBeIsolated[i];
0029
0030 if (verbose_)
0031 cout << "to be isolated: " << toBeIso << endl;
0032
0033 isoDeposits_.push_back(buildIsoDeposit(toBeIso, forIsolation));
0034 }
0035
0036 if (verbose_) {
0037 cout << "PFIsoDepositAlgo " << endl;
0038 }
0039
0040 return isoDeposits_;
0041 }
0042
0043 IsoDeposit PFIsoDepositAlgo::buildIsoDeposit(const Particle& particle, const ParticleCollection& forIsolation) const {
0044 reco::isodeposit::Direction pfDir(particle.eta(), particle.phi());
0045
0046
0047
0048
0049 IsoDeposit isoDep(pfDir);
0050
0051 for (unsigned i = 0; i < forIsolation.size(); i++) {
0052 const reco::PFCandidate& pfc = forIsolation[i];
0053
0054
0055
0056 if (sameParticle(particle, pfc))
0057 continue;
0058
0059 const XYZTLorentzVector& pvi(pfc.p4());
0060 reco::isodeposit::Direction dirPfc(pfc.eta(), pfc.phi());
0061 double dR = pfDir.deltaR(dirPfc);
0062
0063
0064 double maxDeltaRForIsoDep_ = 1;
0065 if (dR > maxDeltaRForIsoDep_) {
0066
0067 continue;
0068 }
0069
0070
0071 if (verbose_)
0072 cout << "\t" << pfc << endl;
0073
0074 double pt = pvi.Pt();
0075 isoDep.addDeposit(dirPfc, pt);
0076 }
0077
0078 return isoDep;
0079 }
0080
0081 bool PFIsoDepositAlgo::sameParticle(const Particle& particle1, const Particle& particle2) const {
0082 double smallNumber = 1e-15;
0083
0084 if (particle1.particleId() != particle2.particleId())
0085 return false;
0086 else if (fabs(particle1.energy() - particle2.energy()) > smallNumber)
0087 return false;
0088 else if (fabs(particle1.eta() - particle2.eta()) > smallNumber)
0089 return false;
0090 else if (fabs(particle1.eta() - particle2.eta()) > smallNumber)
0091 return false;
0092 else
0093 return true;
0094 }