CMS 3D CMS Logo

ShiftedPFCandidateProducerByMatchedObject.cc
Go to the documentation of this file.
2 
3 const double dRDefault = 1000;
4 
6 {
7  srcPFCandidates_ = consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("srcPFCandidates"));
8  srcUnshiftedObjects_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("srcUnshiftedObjects"));
9  srcShiftedObjects_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("srcShiftedObjects"));
10 
11  dRmatch_PFCandidate_ = cfg.getParameter<double>("dRmatch_PFCandidate");
13  dRmatch_Object_ = cfg.exists("dRmatch_Object") ?
14  cfg.getParameter<double>("dRmatch_Object") : 0.1;
16  produces<reco::PFCandidateCollection>();
17 }
18 
20 {
21 // nothing to be done yet...
22 }
23 
25 {
26  edm::Handle<reco::PFCandidateCollection> originalPFCandidates;
27  evt.getByToken(srcPFCandidates_, originalPFCandidates);
28 
30 
31  edm::Handle<CandidateView> unshiftedObjects;
32  evt.getByToken(srcUnshiftedObjects_, unshiftedObjects);
33 
34  edm::Handle<CandidateView> shiftedObjects;
35  evt.getByToken(srcShiftedObjects_, shiftedObjects);
36 
37  objects_.clear();
38 
39  CandidateView::const_iterator shiftedObjectP4_matched;
40  bool isMatched_Object = false;
41  double dR2bestMatch_Object = dRDefault;
42  for ( CandidateView::const_iterator unshiftedObject = unshiftedObjects->begin();
43  unshiftedObject != unshiftedObjects->end(); ++unshiftedObject ) {
44  isMatched_Object = false;
45  dR2bestMatch_Object = dRDefault;
46 
47  for ( CandidateView::const_iterator shiftedObject = shiftedObjects->begin();
48  shiftedObject != shiftedObjects->end(); ++shiftedObject ) {
49  double dR2 = deltaR2(unshiftedObject->p4(), shiftedObject->p4());
50  if ( dR2 < dR2match_Object_ && dR2 < dR2bestMatch_Object ) {
51  shiftedObjectP4_matched = shiftedObject;
52  isMatched_Object = true;
53  dR2bestMatch_Object = dR2;
54  }
55  }
56  if ( isMatched_Object ) {
57  objects_.push_back(objectEntryType(shiftedObjectP4_matched->p4(), unshiftedObject->p4(), sqrt(dR2bestMatch_Object) ));
58  }
59  }
60 
61  auto shiftedPFCandidates = std::make_unique<reco::PFCandidateCollection>();
62 
63  for ( reco::PFCandidateCollection::const_iterator originalPFCandidate = originalPFCandidates->begin();
64  originalPFCandidate != originalPFCandidates->end(); ++originalPFCandidate ) {
65 
66  double shift = 0.;
67  bool applyShift = false;
68  double dR2bestMatch_PFCandidate = dRDefault;
69  for ( std::vector<objectEntryType>::const_iterator object = objects_.begin();
70  object != objects_.end(); ++object ) {
71  if ( !object->isValidMatch_ ) continue;
72  double dR2 = deltaR2(originalPFCandidate->p4(), object->unshiftedObjectP4_);
73  if ( dR2 < dR2match_PFCandidate_ && dR2 < dR2bestMatch_PFCandidate ) {
74  shift = object->shift_;
75  applyShift = true;
76  dR2bestMatch_PFCandidate = dR2;
77  }
78  }
79 
80  reco::Candidate::LorentzVector shiftedPFCandidateP4 = originalPFCandidate->p4();
81  if ( applyShift ) {
82  double shiftedPx = (1. + shift)*originalPFCandidate->px();
83  double shiftedPy = (1. + shift)*originalPFCandidate->py();
84  double shiftedPz = (1. + shift)*originalPFCandidate->pz();
85  double mass = originalPFCandidate->mass();
86  double shiftedEn = sqrt(shiftedPx*shiftedPx + shiftedPy*shiftedPy + shiftedPz*shiftedPz + mass*mass);
87  shiftedPFCandidateP4.SetPxPyPzE(shiftedPx, shiftedPy, shiftedPz, shiftedEn);
88  }
89 
90  reco::PFCandidate shiftedPFCandidate(*originalPFCandidate);
91  shiftedPFCandidate.setP4(shiftedPFCandidateP4);
92 
93  shiftedPFCandidates->push_back(shiftedPFCandidate);
94  }
95 
96  evt.put(std::move(shiftedPFCandidates));
97 }
98 
100 
102 
103 
104 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< edm::View< reco::Candidate > > srcShiftedObjects_
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< edm::View< reco::Candidate > > srcUnshiftedObjects_
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCandidates_
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
static unsigned int const shift
def move(src, dest)
Definition: eostools.py:510
edm::View< Candidate > CandidateView
view of a collection containing candidates
Definition: CandidateFwd.h:23