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