CMS 3D CMS Logo

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