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 // $Id: HiGenCleaner.cc,v 1.3 2010/06/01 14:44:22 yilmaz Exp $
17 //
18 //
19 
20 // 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 template <class T2>
52 class HiGenCleaner : public edm::EDProducer {
53 public:
54  typedef std::vector<T2> T2Collection;
55  explicit HiGenCleaner(const edm::ParameterSet&);
56  ~HiGenCleaner();
57 
58  private:
59  virtual void produce(edm::Event&, const edm::EventSetup&);
60  // ----------member data ---------------------------
61 
63  double deltaR_;
64  double ptCut_;
65  bool makeNew_;
66  bool fillDummy_;
67 
68 };
69 
70 //
71 // constants, enums and typedefs
72 //
73 
74 
75 //
76 // static data member definitions
77 //
78 
79 //
80 // constructors and destructor
81 //
82 
83 template <class T2>
85  jetSrc_(iConfig.getParameter<InputTag>( "src")),
86  deltaR_(iConfig.getParameter<double>("deltaR")),
87  ptCut_(iConfig.getParameter<double>("ptCut")),
88  makeNew_(iConfig.getUntrackedParameter<bool>("createNewCollection",true)),
89  fillDummy_(iConfig.getUntrackedParameter<bool>("fillDummyEntries",true))
90 {
91  std::string alias = jetSrc_.label();
92  produces<T2Collection>().setBranchAlias (alias);
93 }
94 
95 template <class T2>
97 {
98  // do anything here that needs to be done at desctruction time
99  // (e.g. close files, deallocate resources etc.)
100 }
101 
102 
103 //
104 // member functions
105 //
106 
107 // ------------ method called to produce the data ------------
108 template <class T2>
109 void
111 {
112  using namespace edm;
113  using namespace reco;
114 
115  auto_ptr<T2Collection> jets;
116  jets = auto_ptr<T2Collection>(new T2Collection);
117 
118  edm::Handle<edm::View<T2> > genjets;
119  iEvent.getByLabel(jetSrc_,genjets);
120 
121  int jetsize = genjets->size();
122 
123  vector<int> selection;
124  for(int ijet = 0; ijet < jetsize; ++ijet){
125  selection.push_back(-1);
126  }
127 
128  vector<int> selectedIndices;
129  vector<int> removedIndices;
130 
131  for(int ijet = 0; ijet < jetsize; ++ijet){
132 
133  const T2* jet1 = &((*genjets)[ijet]);
134 
135  if(selection[ijet] == -1){
136  selection[ijet] = 1;
137  for(int ijet2 = 0; ijet2 < jetsize; ++ijet2){
138 
139  if(ijet2 == ijet) continue;
140 
141  const T2* jet2 = &((*genjets)[ijet2]);
142 
143  if(Geom::deltaR(jet1->momentum(),jet2->momentum()) < deltaR_){
144  if(jet1->et() < jet2->et()){
145  selection[ijet] = 0;
146  removedIndices.push_back(ijet);
147  break;
148  }else{
149  selection[ijet2] = 0;
150  removedIndices.push_back(ijet2);
151  }
152  }
153  }
154  }
155 
156  double etjet = ((*genjets)[ijet]).et();
157 
158  if(selection[ijet] == 1 && etjet > ptCut_){
159  selectedIndices.push_back(ijet);
160  jets->push_back(*jet1);
161  }
162  }
163  iEvent.put(jets);
164 
165 }
166 
169 
172 
selection
main part
Definition: corrVsCorr.py:98
double deltaR_
std::vector< T2 > T2Collection
Definition: HiGenCleaner.cc:54
DEFINE_FWK_MODULE(HiMixingModule)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
vector< PseudoJet > jets
HiGenCleaner< reco::GenParticle > HiPartonCleaner
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
InputTag jetSrc_
Definition: HiGenCleaner.cc:62
virtual void produce(edm::Event &, const edm::EventSetup &)
HiGenCleaner< reco::GenJet > HiGenJetCleaner
std::string const & label() const
Definition: InputTag.h:25
double deltaR(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:84
HiGenCleaner(const edm::ParameterSet &)
Definition: HiGenCleaner.cc:84
double deltaR_
Definition: HiGenCleaner.cc:63