CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ShiftedJetProducerByMatchedObject.cc
Go to the documentation of this file.
1 
22 
23 #include <string>
24 #include <vector>
25 
26 template <typename T>
28  typedef std::vector<T> JetCollection;
29 
30 public:
33 
34 private:
35  void produce(edm::Event&, const edm::EventSetup&) override;
36 
38 
42 
43  double dRmatch_Jet_;
45 
46  double dR2match_Jet_;
48 
49  struct objectEntryType {
51  const reco::Candidate::LorentzVector& unshiftedObjectP4,
52  double dRmatch)
53  : shiftedObjectP4_(shiftedObjectP4),
54  unshiftedObjectP4_(unshiftedObjectP4),
55  dRmatch_(dRmatch),
57  if (unshiftedObjectP4.energy() > 0.) {
58  shift_ = (shiftedObjectP4.energy() / unshiftedObjectP4.energy()) - 1.;
59  isValidMatch_ = true;
60  }
61  }
65  double shift_;
66  double dRmatch_;
68  };
69 
70  std::vector<objectEntryType> objects_;
71 };
72 
73 template <typename T>
75  : moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
76  srcJets_ = consumes<JetCollection>(cfg.getParameter<edm::InputTag>("srcJets"));
77  srcUnshiftedObjects_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("srcUnshiftedObjects"));
78  srcShiftedObjects_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("srcShiftedObjects"));
79 
80  dRmatch_Jet_ = cfg.getParameter<double>("dRmatch_Jet");
81  dRmatch_Object_ = cfg.exists("dRmatch_Object") ? cfg.getParameter<double>("dRmatch_Object") : 0.1;
82 
84  dR2match_Object_ = dRmatch_Object_ * dRmatch_Object_;
85 
86  produces<JetCollection>();
87 }
88 
89 template <typename T>
91  // nothing to be done yet...
92 }
93 
94 template <typename T>
96  edm::Handle<JetCollection> originalJets;
97  evt.getByToken(srcJets_, originalJets);
98 
99  edm::Handle<reco::CandidateView> unshiftedObjects;
100  evt.getByToken(srcUnshiftedObjects_, unshiftedObjects);
101 
102  edm::Handle<reco::CandidateView> shiftedObjects;
103  evt.getByToken(srcShiftedObjects_, shiftedObjects);
104 
105  objects_.clear();
106 
107  std::vector<bool> match(shiftedObjects->size(), false);
108  int prevMatch = -1;
109  int cnt = 0;
110 
111  for (reco::CandidateView::const_iterator unshiftedObject = unshiftedObjects->begin();
112  unshiftedObject != unshiftedObjects->end();
113  ++unshiftedObject) {
114  bool isMatched_Object = false;
115  double dR2bestMatch_Object = std::numeric_limits<double>::max();
116  prevMatch = -1;
117  cnt = 0;
118 
119  reco::Candidate::LorentzVector shiftedObjectP4_matched;
120  for (reco::CandidateView::const_iterator shiftedObject = shiftedObjects->begin();
121  shiftedObject != shiftedObjects->end();
122  ++shiftedObject) {
123  if (match[cnt])
124  continue;
125 
126  double dR2 = deltaR2(unshiftedObject->p4(), shiftedObject->p4());
127  if (dR2 < dR2match_Object_ && dR2 < dR2bestMatch_Object) {
128  shiftedObjectP4_matched = shiftedObject->p4();
129  isMatched_Object = true;
130  dR2bestMatch_Object = dR2;
131 
132  prevMatch = cnt;
133  }
134  cnt++;
135  }
136  if (isMatched_Object) {
137  //Ambiguity removal
138  match[prevMatch] = true;
139  objects_.push_back(objectEntryType(shiftedObjectP4_matched, unshiftedObject->p4(), sqrt(dR2bestMatch_Object)));
140  }
141  }
142 
143  match.assign(objects_.size(), false);
144 
145  auto shiftedJets = std::make_unique<JetCollection>();
146 
147  for (typename JetCollection::const_iterator originalJet = originalJets->begin(); originalJet != originalJets->end();
148  ++originalJet) {
149  double shift = 0.;
150  bool applyShift = false;
151  double dR2bestMatch_Jet = std::numeric_limits<double>::max();
152  prevMatch = -1;
153  cnt = 0;
154 
155  for (typename std::vector<objectEntryType>::const_iterator object = objects_.begin(); object != objects_.end();
156  ++object) {
157  if (!object->isValidMatch_)
158  continue;
159  if (match[cnt])
160  continue;
161 
162  double dR2 = deltaR2(originalJet->p4(), object->unshiftedObjectP4_);
163  if (dR2 < dR2match_Jet_ && dR2 < dR2bestMatch_Jet) {
164  shift = object->shift_;
165  applyShift = true;
166  dR2bestMatch_Jet = dR2;
167 
168  prevMatch = cnt;
169  }
170  cnt++;
171  }
172 
173  reco::Candidate::LorentzVector shiftedJetP4 = originalJet->p4();
174  if (applyShift) {
175  //Ambiguity removal
176  match[prevMatch] = true;
177 
178  shiftedJetP4 *= (1. + shift);
179  }
180 
181  T shiftedJet(*originalJet);
182  shiftedJet.setP4(shiftedJetP4);
183 
184  shiftedJets->push_back(shiftedJet);
185  }
186 
187  evt.put(std::move(shiftedJets));
188 }
189 
192 
195 
197 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple cfg
Definition: looper.py:296
ShiftedJetProducerByMatchedObjectT< reco::CaloJet > ShiftedCaloJetProducerByMatchedObject
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void produce(edm::Event &, const edm::EventSetup &) override
bool exists(std::string const &parameterName) const
checks if a parameter exists
T sqrt(T t)
Definition: SSEVec.h:19
def move
Definition: eostools.py:511
edm::EDGetTokenT< edm::View< reco::Candidate > > srcShiftedObjects_
ShiftedJetProducerByMatchedObjectT< reco::PFJet > ShiftedPFJetProducerByMatchedObject
edm::EDGetTokenT< edm::View< reco::Candidate > > srcUnshiftedObjects_
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:36
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
static unsigned int const shift
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
moduleLabel_(iConfig.getParameter< string >("@module_label"))
long double T
ShiftedJetProducerByMatchedObjectT(const edm::ParameterSet &)
objectEntryType(const reco::Candidate::LorentzVector &shiftedObjectP4, const reco::Candidate::LorentzVector &unshiftedObjectP4, double dRmatch)