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
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 jets = std::make_unique<T2Collection>();
97 
98  edm::Handle<edm::View<T2> > genjets;
99  iEvent.getByToken(jetSrc_,genjets);
100 
101  int jetsize = genjets->size();
102 
103  vector<int> selection;
104  for(int ijet = 0; ijet < jetsize; ++ijet){
105  selection.push_back(-1);
106  }
107 
108  vector<int> selectedIndices;
109  vector<int> removedIndices;
110 
111  for(int ijet = 0; ijet < jetsize; ++ijet){
112 
113  const T2* jet1 = &((*genjets)[ijet]);
114 
115  if(selection[ijet] == -1){
116  selection[ijet] = 1;
117  for(int ijet2 = 0; ijet2 < jetsize; ++ijet2){
118 
119  if(ijet2 == ijet) continue;
120 
121  const T2* jet2 = &((*genjets)[ijet2]);
122 
123  if(Geom::deltaR(jet1->momentum(),jet2->momentum()) < deltaR_){
124  if(jet1->et() < jet2->et()){
125  selection[ijet] = 0;
126  removedIndices.push_back(ijet);
127  break;
128  }else{
129  selection[ijet2] = 0;
130  removedIndices.push_back(ijet2);
131  }
132  }
133  }
134  }
135 
136  double etjet = ((*genjets)[ijet]).et();
137 
138  if(selection[ijet] == 1 && etjet > ptCut_){
139  selectedIndices.push_back(ijet);
140  jets->push_back(*jet1);
141  }
142  }
143  iEvent.put(std::move(jets));
144 
145 }
146 
147 
148 template class HiGenCleaner<reco::GenParticle>;
149 template class HiGenCleaner<reco::GenJet>;
double ptCut_
Definition: HiGenCleaner.h:53
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
selection
main part
Definition: corrVsCorr.py:98
double deltaR_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< edm::View< T2 > > jetSrc_
Definition: HiGenCleaner.h:51
vector< PseudoJet > jets
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
et
define resolution functions of each parameter
fixed size matrix
HLT enums.
HiGenCleaner(const edm::ParameterSet &)
Definition: HiGenCleaner.cc:65
def move(src, dest)
Definition: eostools.py:510
~HiGenCleaner() override
Definition: HiGenCleaner.cc:77
double deltaR_
Definition: HiGenCleaner.h:52