CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ShiftedPFCandidateProducerForPFMEtMVA.cc
Go to the documentation of this file.
2 
4 
6  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
7  , srcPFCandidatesToken_(consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("srcPFCandidates")))
8  , srcUnshiftedObjectsToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcUnshiftedObjects")))
9  , srcShiftedObjectsToken_(consumes<CandidateView>(cfg.getParameter<edm::InputTag>("srcShiftedObjects")))
10 {
11  dRmatch_PFCandidate_ = cfg.getParameter<double>("dRmatch_PFCandidate");
12  dRmatch_Object_ = cfg.exists("dRmatch_Object") ?
13  cfg.getParameter<double>("dRmatch_Object") : 0.1;
14 
15  produces<reco::PFCandidateCollection>();
16 }
17 
19 {
20 // nothing to be done yet...
21 }
22 
24 {
25  edm::Handle<reco::PFCandidateCollection> originalPFCandidates;
26  evt.getByToken(srcPFCandidatesToken_, originalPFCandidates);
27 
28  edm::Handle<CandidateView> unshiftedObjects;
29  evt.getByToken(srcUnshiftedObjectsToken_, unshiftedObjects);
30 
31  edm::Handle<CandidateView> shiftedObjects;
32  evt.getByToken(srcShiftedObjectsToken_, shiftedObjects);
33 
34  objects_.clear();
35 
36  for ( CandidateView::const_iterator unshiftedObject = unshiftedObjects->begin();
37  unshiftedObject != unshiftedObjects->end(); ++unshiftedObject ) {
38  bool isMatched_Object = false;
39  double dRbestMatch_Object = 1.e+3;
40  reco::Candidate::LorentzVector shiftedObjectP4_matched;
41  for ( CandidateView::const_iterator shiftedObject = shiftedObjects->begin();
42  shiftedObject != shiftedObjects->end(); ++shiftedObject ) {
43  double dR = deltaR(unshiftedObject->p4(), shiftedObject->p4());
44  if ( dR < dRmatch_Object_ && dR < dRbestMatch_Object ) {
45  shiftedObjectP4_matched = shiftedObject->p4();
46  isMatched_Object = true;
47  dRbestMatch_Object = dR;
48  }
49  }
50  if ( isMatched_Object ) {
51  objects_.push_back(objectEntryType(shiftedObjectP4_matched, unshiftedObject->p4(), dRbestMatch_Object));
52  }
53  }
54 
55  std::auto_ptr<reco::PFCandidateCollection> shiftedPFCandidates(new reco::PFCandidateCollection);
56 
57  for ( reco::PFCandidateCollection::const_iterator originalPFCandidate = originalPFCandidates->begin();
58  originalPFCandidate != originalPFCandidates->end(); ++originalPFCandidate ) {
59 
60  double shift = 0.;
61  bool applyShift = false;
62  double dRbestMatch_PFCandidate = 1.e+3;
63  for ( std::vector<objectEntryType>::const_iterator object = objects_.begin();
64  object != objects_.end(); ++object ) {
65  if ( !object->isValidMatch_ ) continue;
66  double dR = deltaR(originalPFCandidate->p4(), object->unshiftedObjectP4_);
67  if ( dR < dRmatch_PFCandidate_ && dR < dRbestMatch_PFCandidate ) {
68  shift = object->shift_;
69  applyShift = true;
70  dRbestMatch_PFCandidate = dR;
71  }
72  }
73 
74  reco::Candidate::LorentzVector shiftedPFCandidateP4 = originalPFCandidate->p4();
75  if ( applyShift ) shiftedPFCandidateP4 *= (1. + shift);
76 
77  reco::PFCandidate shiftedPFCandidate(*originalPFCandidate);
78  shiftedPFCandidate.setP4(shiftedPFCandidateP4);
79 
80  shiftedPFCandidates->push_back(shiftedPFCandidate);
81  }
82 
83  evt.put(shiftedPFCandidates);
84 }
85 
87 
89 
90 
91 
T getParameter(std::string const &) const
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
edm::EDGetTokenT< CandidateView > srcUnshiftedObjectsToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual void setP4(const LorentzVector &p4)
set 4-momentum
void produce(edm::Event &, const edm::EventSetup &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
std::vector< reco::PFCandidate > PFCandidateCollection
collection of PFCandidates
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
list object
Definition: dbtoconf.py:77
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:38
static unsigned int const shift
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCandidatesToken_
moduleLabel_(iConfig.getParameter< string >("@module_label"))