CMS 3D CMS Logo

ShiftedPFCandidateProducerForPFMVAMEt.cc
Go to the documentation of this file.
2 
4 
5 const double dRDefault = 1000;
6 
8  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
9  , srcPFCandidatesToken_(consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("srcPFCandidates")))
10  , srcUnshiftedObjectsToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcUnshiftedObjects")))
11  , srcShiftedObjectsToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcShiftedObjects")))
12 {
13  dRmatch_PFCandidate_ = cfg.getParameter<double>("dRmatch_PFCandidate");
14  dRmatch_Object_ = cfg.exists("dRmatch_Object") ?
15  cfg.getParameter<double>("dRmatch_Object") : 0.1;
16 
18  dR2match_Object_ = dRmatch_Object_*dRmatch_Object_;
19 
20  produces<reco::PFCandidateCollection>();
21 }
22 
24 {
25 // nothing to be done yet...
26 }
27 
29 {
30  edm::Handle<reco::PFCandidateCollection> originalPFCandidates;
31  evt.getByToken(srcPFCandidatesToken_, originalPFCandidates);
32 
33  edm::Handle<CandidateView> unshiftedObjects;
34  evt.getByToken(srcUnshiftedObjectsToken_, unshiftedObjects);
35 
36  edm::Handle<CandidateView> shiftedObjects;
37  evt.getByToken(srcShiftedObjectsToken_, shiftedObjects);
38 
39  objects_.clear();
40 
41  std::vector<bool> match(shiftedObjects->size(), false);
42  int prevMatch=-1;
43  int cnt = 0;
44 
45  for ( CandidateView::const_iterator unshiftedObject = unshiftedObjects->begin();
46  unshiftedObject != unshiftedObjects->end(); ++unshiftedObject ) {
47  bool isMatched_Object = false;
48  double dR2bestMatch_Object = dRDefault;
49  reco::Candidate::LorentzVector shiftedObjectP4_matched;
50  prevMatch=-1;
51 
52  for ( CandidateView::const_iterator shiftedObject = shiftedObjects->begin();
53  shiftedObject != shiftedObjects->end(); ++shiftedObject ) {
54  if( match[ cnt ] ) continue;
55 
56  double dR2 = deltaR2(unshiftedObject->p4(), shiftedObject->p4());
57  if ( dR2 < dR2match_Object_ && dR2 < dR2bestMatch_Object ) {
58  shiftedObjectP4_matched = shiftedObject->p4();
59  isMatched_Object = true;
60  dR2bestMatch_Object = dR2;
61 
62  prevMatch = cnt;
63  }
64  cnt++;
65  }
66  if ( isMatched_Object ) {
67  //Ambiguity removal
68  match[ prevMatch ] = true;
69  objects_.push_back(objectEntryType(shiftedObjectP4_matched, unshiftedObject->p4(), sqrt(dR2bestMatch_Object)));
70  }
71  }
72 
73  match.assign(objects_.size(), false);
74 
75  auto shiftedPFCandidates = std::make_unique<reco::PFCandidateCollection>();
76 
77  for ( reco::PFCandidateCollection::const_iterator originalPFCandidate = originalPFCandidates->begin();
78  originalPFCandidate != originalPFCandidates->end(); ++originalPFCandidate ) {
79 
80  double shift = 0.;
81  bool applyShift = false;
82  double dR2bestMatch_PFCandidate = dRDefault;
83  prevMatch=-1;
84  cnt = 0;
85 
86  for ( std::vector<objectEntryType>::const_iterator object = objects_.begin();
87  object != objects_.end(); ++object ) {
88  if ( !object->isValidMatch_ ) continue;
89  if( match[ cnt ] ) continue;
90 
91  double dR2 = deltaR2(originalPFCandidate->p4(), object->unshiftedObjectP4_);
92  if ( dR2 < dR2match_PFCandidate_ && dR2 < dR2bestMatch_PFCandidate ) {
93  shift = object->shift_;
94  applyShift = true;
95  dR2bestMatch_PFCandidate = dR2;
96 
97  prevMatch = cnt;
98  }
99  cnt++;
100  }
101 
102  reco::Candidate::LorentzVector shiftedPFCandidateP4 = originalPFCandidate->p4();
103  if ( applyShift ) {
104  //Ambiguity removal
105  match[ prevMatch ] = true;
106 
107  shiftedPFCandidateP4 *= (1. + shift);
108  }
109 
110  reco::PFCandidate shiftedPFCandidate(*originalPFCandidate);
111  shiftedPFCandidate.setP4(shiftedPFCandidateP4);
112 
113  shiftedPFCandidates->push_back(shiftedPFCandidate);
114  }
115 
116  evt.put(std::move(shiftedPFCandidates));
117 }
118 
120 
122 
123 
124 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
void produce(edm::Event &, const edm::EventSetup &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
size_type size() const
const double dRDefault
const_iterator begin() const
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
fixed size matrix
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
static unsigned int const shift
const_iterator end() const
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< CandidateView > srcUnshiftedObjectsToken_
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCandidatesToken_