CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
27 
30 
33 
39 
42 
43 using namespace std;
44 using namespace edm;
45 
46 //
47 // class decleration
48 //
49 
50 //
51 // constants, enums and typedefs
52 //
53 
54 //
55 // static data member definitions
56 //
57 
58 //
59 // constructors and destructor
60 //
61 
62 template <class T2>
64  : jetSrc_(consumes<edm::View<T2> >(iConfig.getParameter<edm::InputTag>("src"))),
65  deltaR_(iConfig.getParameter<double>("deltaR")),
66  ptCut_(iConfig.getParameter<double>("ptCut")),
67  makeNew_(iConfig.getUntrackedParameter<bool>("createNewCollection", true)),
68  fillDummy_(iConfig.getUntrackedParameter<bool>("fillDummyEntries", true)) {
69  std::string alias = (iConfig.getParameter<InputTag>("src")).label();
70  produces<T2Collection>().setBranchAlias(alias);
71 }
72 
73 template <class T2>
75  // do anything here that needs to be done at desctruction time
76  // (e.g. close files, deallocate resources etc.)
77 }
78 
79 //
80 // member functions
81 //
82 
83 // ------------ method called to produce the data ------------
84 template <class T2>
86  using namespace edm;
87  using namespace reco;
88 
89  auto jets = std::make_unique<T2Collection>();
90 
91  edm::Handle<edm::View<T2> > genjets;
92  iEvent.getByToken(jetSrc_, genjets);
93 
94  int jetsize = genjets->size();
95 
96  vector<int> selection;
97  selection.reserve(jetsize);
98  for (int ijet = 0; ijet < jetsize; ++ijet) {
99  selection.push_back(-1);
100  }
101 
102  vector<int> selectedIndices;
103  vector<int> removedIndices;
104 
105  for (int ijet = 0; ijet < jetsize; ++ijet) {
106  const T2* jet1 = &((*genjets)[ijet]);
107 
108  if (selection[ijet] == -1) {
109  selection[ijet] = 1;
110  for (int ijet2 = 0; ijet2 < jetsize; ++ijet2) {
111  if (ijet2 == ijet)
112  continue;
113 
114  const T2* jet2 = &((*genjets)[ijet2]);
115 
116  if (Geom::deltaR(jet1->momentum(), jet2->momentum()) < deltaR_) {
117  if (jet1->et() < jet2->et()) {
118  selection[ijet] = 0;
119  removedIndices.push_back(ijet);
120  break;
121  } else {
122  selection[ijet2] = 0;
123  removedIndices.push_back(ijet2);
124  }
125  }
126  }
127  }
128 
129  double etjet = ((*genjets)[ijet]).et();
130 
131  if (selection[ijet] == 1 && etjet > ptCut_) {
132  selectedIndices.push_back(ijet);
133  jets->push_back(*jet1);
134  }
135  }
136  iEvent.put(std::move(jets));
137 }
138 
139 template class HiGenCleaner<reco::GenParticle>;
140 template class HiGenCleaner<reco::GenJet>;
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
selection
main part
Definition: corrVsCorr.py:100
int iEvent
Definition: GenABIO.cc:224
vector< PseudoJet > jets
def move
Definition: eostools.py:511
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: HiGenCleaner.cc:85
HiGenCleaner(const edm::ParameterSet &)
Definition: HiGenCleaner.cc:63
~HiGenCleaner() override
Definition: HiGenCleaner.cc:74