CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ShiftedPFCandidateProducerForPFMVAMEt.cc
Go to the documentation of this file.
1 
26 
27 #include <string>
28 #include <vector>
29 
31 public:
34 
35 private:
37 
38  void produce(edm::Event&, const edm::EventSetup&) override;
39 
41 
45 
48 
51 
52  struct objectEntryType {
54  const reco::Candidate::LorentzVector& unshiftedObjectP4,
55  double dRmatch)
56  : shiftedObjectP4_(shiftedObjectP4),
57  unshiftedObjectP4_(unshiftedObjectP4),
58  dRmatch_(dRmatch),
60  if (unshiftedObjectP4.energy() > 0.) {
61  shift_ = (shiftedObjectP4.energy() / unshiftedObjectP4.energy()) - 1.;
62  isValidMatch_ = true;
63  }
64  }
68  double dRmatch_;
69  double shift_;
71  };
72 
73  std::vector<objectEntryType> objects_;
74 };
75 
76 const double dRDefault = 1000;
77 
79  : moduleLabel_(cfg.getParameter<std::string>("@module_label")),
80  srcPFCandidatesToken_(consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("srcPFCandidates"))),
81  srcUnshiftedObjectsToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcUnshiftedObjects"))),
82  srcShiftedObjectsToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcShiftedObjects"))) {
83  dRmatch_PFCandidate_ = cfg.getParameter<double>("dRmatch_PFCandidate");
84  dRmatch_Object_ = cfg.exists("dRmatch_Object") ? cfg.getParameter<double>("dRmatch_Object") : 0.1;
85 
87  dR2match_Object_ = dRmatch_Object_ * dRmatch_Object_;
88 
89  produces<reco::PFCandidateCollection>();
90 }
91 
93  // nothing to be done yet...
94 }
95 
97  edm::Handle<reco::PFCandidateCollection> originalPFCandidates;
98  evt.getByToken(srcPFCandidatesToken_, originalPFCandidates);
99 
100  edm::Handle<CandidateView> unshiftedObjects;
101  evt.getByToken(srcUnshiftedObjectsToken_, unshiftedObjects);
102 
103  edm::Handle<CandidateView> shiftedObjects;
104  evt.getByToken(srcShiftedObjectsToken_, shiftedObjects);
105 
106  objects_.clear();
107 
108  std::vector<bool> match(shiftedObjects->size(), false);
109  int prevMatch = -1;
110  int cnt = 0;
111 
112  for (CandidateView::const_iterator unshiftedObject = unshiftedObjects->begin();
113  unshiftedObject != unshiftedObjects->end();
114  ++unshiftedObject) {
115  bool isMatched_Object = false;
116  double dR2bestMatch_Object = dRDefault;
117  reco::Candidate::LorentzVector shiftedObjectP4_matched;
118  prevMatch = -1;
119 
120  for (CandidateView::const_iterator shiftedObject = shiftedObjects->begin(); shiftedObject != shiftedObjects->end();
121  ++shiftedObject) {
122  if (match[cnt])
123  continue;
124 
125  double dR2 = deltaR2(unshiftedObject->p4(), shiftedObject->p4());
126  if (dR2 < dR2match_Object_ && dR2 < dR2bestMatch_Object) {
127  shiftedObjectP4_matched = shiftedObject->p4();
128  isMatched_Object = true;
129  dR2bestMatch_Object = dR2;
130 
131  prevMatch = cnt;
132  }
133  cnt++;
134  }
135  if (isMatched_Object) {
136  //Ambiguity removal
137  match[prevMatch] = true;
138  objects_.push_back(objectEntryType(shiftedObjectP4_matched, unshiftedObject->p4(), sqrt(dR2bestMatch_Object)));
139  }
140  }
141 
142  match.assign(objects_.size(), false);
143 
144  auto shiftedPFCandidates = std::make_unique<reco::PFCandidateCollection>();
145 
146  for (reco::PFCandidateCollection::const_iterator originalPFCandidate = originalPFCandidates->begin();
147  originalPFCandidate != originalPFCandidates->end();
148  ++originalPFCandidate) {
149  double shift = 0.;
150  bool applyShift = false;
151  double dR2bestMatch_PFCandidate = dRDefault;
152  prevMatch = -1;
153  cnt = 0;
154 
155  for (std::vector<objectEntryType>::const_iterator object = objects_.begin(); object != objects_.end(); ++object) {
156  if (!object->isValidMatch_)
157  continue;
158  if (match[cnt])
159  continue;
160 
161  double dR2 = deltaR2(originalPFCandidate->p4(), object->unshiftedObjectP4_);
162  if (dR2 < dR2match_PFCandidate_ && dR2 < dR2bestMatch_PFCandidate) {
163  shift = object->shift_;
164  applyShift = true;
165  dR2bestMatch_PFCandidate = dR2;
166 
167  prevMatch = cnt;
168  }
169  cnt++;
170  }
171 
172  reco::Candidate::LorentzVector shiftedPFCandidateP4 = originalPFCandidate->p4();
173  if (applyShift) {
174  //Ambiguity removal
175  match[prevMatch] = true;
176 
177  shiftedPFCandidateP4 *= (1. + shift);
178  }
179 
180  reco::PFCandidate shiftedPFCandidate(*originalPFCandidate);
181  shiftedPFCandidate.setP4(shiftedPFCandidateP4);
182 
183  shiftedPFCandidates->push_back(shiftedPFCandidate);
184  }
185 
186  evt.put(std::move(shiftedPFCandidates));
187 }
188 
190 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple cfg
Definition: looper.py:296
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< l1t::PFCandidate > PFCandidateCollection
Definition: PFCandidate.h:57
bool exists(std::string const &parameterName) const
checks if a parameter exists
objectEntryType(const reco::Candidate::LorentzVector &shiftedObjectP4, const reco::Candidate::LorentzVector &unshiftedObjectP4, double dRmatch)
T sqrt(T t)
Definition: SSEVec.h:19
def move
Definition: eostools.py:511
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
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
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
void produce(edm::Event &, const edm::EventSetup &) override
moduleLabel_(iConfig.getParameter< string >("@module_label"))
void setP4(const LorentzVector &p4) final
set 4-momentum
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCandidatesToken_