CMS 3D CMS Logo

HiGenCleaner.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HiGenCleaner
4 // Class: HiGenCleaner
5 //
13 //
14 // Original Author: Yetkin Yilmaz
15 // Created: Tue Jul 21 04:26:01 EDT 2009
16 //
17 //
18 
19 // system include files
21 #include <memory>
22 #include <vector>
23 
24 // user include files
26 
29 
32 
38 
41 
42 using namespace std;
43 using namespace edm;
44 
45 //
46 // class decleration
47 //
48 
49 //
50 // constants, enums and typedefs
51 //
52 
53 //
54 // static data member definitions
55 //
56 
57 //
58 // constructors and destructor
59 //
60 
61 template <class T2>
63  : jetSrc_(consumes<edm::View<T2> >(iConfig.getParameter<edm::InputTag>("src"))),
64  deltaR_(iConfig.getParameter<double>("deltaR")),
65  ptCut_(iConfig.getParameter<double>("ptCut")),
66  makeNew_(iConfig.getUntrackedParameter<bool>("createNewCollection", true)),
67  fillDummy_(iConfig.getUntrackedParameter<bool>("fillDummyEntries", true)) {
68  std::string alias = (iConfig.getParameter<InputTag>("src")).label();
69  produces<T2Collection>().setBranchAlias(alias);
70 }
71 
72 template <class T2>
74  // do anything here that needs to be done at desctruction time
75  // (e.g. close files, deallocate resources etc.)
76 }
77 
78 //
79 // member functions
80 //
81 
82 // ------------ method called to produce the data ------------
83 template <class T2>
85  using namespace edm;
86  using namespace reco;
87 
88  auto jets = std::make_unique<T2Collection>();
89 
90  edm::Handle<edm::View<T2> > genjets;
91  iEvent.getByToken(jetSrc_, genjets);
92 
93  int jetsize = genjets->size();
94 
95  vector<int> selection;
96  selection.reserve(jetsize);
97  for (int ijet = 0; ijet < jetsize; ++ijet) {
98  selection.push_back(-1);
99  }
100 
101  vector<int> selectedIndices;
102  vector<int> removedIndices;
103 
104  for (int ijet = 0; ijet < jetsize; ++ijet) {
105  const T2* jet1 = &((*genjets)[ijet]);
106 
107  if (selection[ijet] == -1) {
108  selection[ijet] = 1;
109  for (int ijet2 = 0; ijet2 < jetsize; ++ijet2) {
110  if (ijet2 == ijet)
111  continue;
112 
113  const T2* jet2 = &((*genjets)[ijet2]);
114 
115  if (Geom::deltaR(jet1->momentum(), jet2->momentum()) < deltaR_) {
116  if (jet1->et() < jet2->et()) {
117  selection[ijet] = 0;
118  removedIndices.push_back(ijet);
119  break;
120  } else {
121  selection[ijet2] = 0;
122  removedIndices.push_back(ijet2);
123  }
124  }
125  }
126  }
127 
128  double etjet = ((*genjets)[ijet]).et();
129 
130  if (selection[ijet] == 1 && etjet > ptCut_) {
131  selectedIndices.push_back(ijet);
132  jets->push_back(*jet1);
133  }
134  }
135  iEvent.put(std::move(jets));
136 }
137 
138 template class HiGenCleaner<reco::GenParticle>;
139 template class HiGenCleaner<reco::GenJet>;
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
selection
main part
Definition: corrVsCorr.py:100
int iEvent
Definition: GenABIO.cc:224
fixed size matrix
HLT enums.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: HiGenCleaner.cc:84
HiGenCleaner(const edm::ParameterSet &)
Definition: HiGenCleaner.cc:62
def move(src, dest)
Definition: eostools.py:511
~HiGenCleaner() override
Definition: HiGenCleaner.cc:73