File indexing completed on 2023-10-25 09:59:36
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include "RecoHI/HiJetAlgos/interface/HiGenCleaner.h"
0021 #include <memory>
0022 #include <vector>
0023
0024
0025 #include "FWCore/Framework/interface/Frameworkfwd.h"
0026
0027 #include "FWCore/Framework/interface/Event.h"
0028 #include "FWCore/Framework/interface/MakerMacros.h"
0029
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 #include "FWCore/Utilities/interface/InputTag.h"
0032
0033 #include "DataFormats/Common/interface/View.h"
0034 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0035 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
0036 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0037 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
0038
0039 #include "DataFormats/Math/interface/Point3D.h"
0040 #include "DataFormats/Math/interface/LorentzVector.h"
0041
0042 using namespace std;
0043 using namespace edm;
0044
0045
0046
0047
0048
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058
0059
0060
0061 template <class T2>
0062 HiGenCleaner<T2>::HiGenCleaner(const edm::ParameterSet& iConfig)
0063 : jetSrc_(consumes<edm::View<T2> >(iConfig.getParameter<edm::InputTag>("src"))),
0064 deltaR_(iConfig.getParameter<double>("deltaR")),
0065 ptCut_(iConfig.getParameter<double>("ptCut")),
0066 makeNew_(iConfig.getUntrackedParameter<bool>("createNewCollection", true)),
0067 fillDummy_(iConfig.getUntrackedParameter<bool>("fillDummyEntries", true)) {
0068 std::string alias = (iConfig.getParameter<InputTag>("src")).label();
0069 produces<T2Collection>().setBranchAlias(alias);
0070 }
0071
0072 template <class T2>
0073 HiGenCleaner<T2>::~HiGenCleaner() {
0074
0075
0076 }
0077
0078
0079
0080
0081
0082
0083 template <class T2>
0084 void HiGenCleaner<T2>::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0085 using namespace edm;
0086 using namespace reco;
0087
0088 auto jets = std::make_unique<T2Collection>();
0089
0090 edm::Handle<edm::View<T2> > genjets;
0091 iEvent.getByToken(jetSrc_, genjets);
0092
0093 int jetsize = genjets->size();
0094
0095 vector<int> selection;
0096 selection.reserve(jetsize);
0097 for (int ijet = 0; ijet < jetsize; ++ijet) {
0098 selection.push_back(-1);
0099 }
0100
0101 vector<int> selectedIndices;
0102 vector<int> removedIndices;
0103
0104 for (int ijet = 0; ijet < jetsize; ++ijet) {
0105 const T2* jet1 = &((*genjets)[ijet]);
0106
0107 if (selection[ijet] == -1) {
0108 selection[ijet] = 1;
0109 for (int ijet2 = 0; ijet2 < jetsize; ++ijet2) {
0110 if (ijet2 == ijet)
0111 continue;
0112
0113 const T2* jet2 = &((*genjets)[ijet2]);
0114
0115 if (Geom::deltaR(jet1->momentum(), jet2->momentum()) < deltaR_) {
0116 if (jet1->et() < jet2->et()) {
0117 selection[ijet] = 0;
0118 removedIndices.push_back(ijet);
0119 break;
0120 } else {
0121 selection[ijet2] = 0;
0122 removedIndices.push_back(ijet2);
0123 }
0124 }
0125 }
0126 }
0127
0128 double etjet = ((*genjets)[ijet]).et();
0129
0130 if (selection[ijet] == 1 && etjet > ptCut_) {
0131 selectedIndices.push_back(ijet);
0132 jets->push_back(*jet1);
0133 }
0134 }
0135 iEvent.put(std::move(jets));
0136 }
0137
0138 template class HiGenCleaner<reco::GenParticle>;
0139 template class HiGenCleaner<reco::GenJet>;