CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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
20 #include <memory>
21 #include <vector>
22 
23 // user include files
26 
29 
32 
38 
41 
42 using namespace std;
43 using namespace edm;
44 
45 
46 //
47 // class decleration
48 //
49 
50 template <class T2>
51 class HiGenCleaner : public edm::EDProducer {
52 public:
53  typedef std::vector<T2> T2Collection;
54  explicit HiGenCleaner(const edm::ParameterSet&);
55  ~HiGenCleaner();
56 
57  private:
58  virtual void produce(edm::Event&, const edm::EventSetup&) override;
59  // ----------member data ---------------------------
60 
62  double deltaR_;
63  double ptCut_;
64  bool makeNew_;
65  bool fillDummy_;
66 
67 };
68 
69 //
70 // constants, enums and typedefs
71 //
72 
73 
74 //
75 // static data member definitions
76 //
77 
78 //
79 // constructors and destructor
80 //
81 
82 template <class T2>
84  jetSrc_(consumes<edm::View<T2> >(iConfig.getParameter<edm::InputTag>("src"))),
85  deltaR_(iConfig.getParameter<double>("deltaR")),
86  ptCut_(iConfig.getParameter<double>("ptCut")),
87  makeNew_(iConfig.getUntrackedParameter<bool>("createNewCollection",true)),
88  fillDummy_(iConfig.getUntrackedParameter<bool>("fillDummyEntries",true))
89 {
90  std::string alias = (iConfig.getParameter<InputTag>( "src")).label();
91  produces<T2Collection>().setBranchAlias (alias);
92 }
93 
94 template <class T2>
96 {
97  // do anything here that needs to be done at desctruction time
98  // (e.g. close files, deallocate resources etc.)
99 }
100 
101 
102 //
103 // member functions
104 //
105 
106 // ------------ method called to produce the data ------------
107 template <class T2>
108 void
110 {
111  using namespace edm;
112  using namespace reco;
113 
114  auto_ptr<T2Collection> jets;
115  jets = auto_ptr<T2Collection>(new T2Collection);
116 
117  edm::Handle<edm::View<T2> > genjets;
118  iEvent.getByToken(jetSrc_,genjets);
119 
120  int jetsize = genjets->size();
121 
122  vector<int> selection;
123  for(int ijet = 0; ijet < jetsize; ++ijet){
124  selection.push_back(-1);
125  }
126 
127  vector<int> selectedIndices;
128  vector<int> removedIndices;
129 
130  for(int ijet = 0; ijet < jetsize; ++ijet){
131 
132  const T2* jet1 = &((*genjets)[ijet]);
133 
134  if(selection[ijet] == -1){
135  selection[ijet] = 1;
136  for(int ijet2 = 0; ijet2 < jetsize; ++ijet2){
137 
138  if(ijet2 == ijet) continue;
139 
140  const T2* jet2 = &((*genjets)[ijet2]);
141 
142  if(Geom::deltaR(jet1->momentum(),jet2->momentum()) < deltaR_){
143  if(jet1->et() < jet2->et()){
144  selection[ijet] = 0;
145  removedIndices.push_back(ijet);
146  break;
147  }else{
148  selection[ijet2] = 0;
149  removedIndices.push_back(ijet2);
150  }
151  }
152  }
153  }
154 
155  double etjet = ((*genjets)[ijet]).et();
156 
157  if(selection[ijet] == 1 && etjet > ptCut_){
158  selectedIndices.push_back(ijet);
159  jets->push_back(*jet1);
160  }
161  }
162  iEvent.put(jets);
163 
164 }
165 
168 
171 
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
selection
main part
Definition: corrVsCorr.py:98
double deltaR_
std::vector< T2 > T2Collection
Definition: HiGenCleaner.cc:53
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::View< T2 > > jetSrc_
Definition: HiGenCleaner.cc:61
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
vector< PseudoJet > jets
HiGenCleaner< reco::GenParticle > HiPartonCleaner
virtual void produce(edm::Event &, const edm::EventSetup &) override
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
HiGenCleaner< reco::GenJet > HiGenJetCleaner
HiGenCleaner(const edm::ParameterSet &)
Definition: HiGenCleaner.cc:83
double deltaR_
Definition: HiGenCleaner.cc:62