CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ShiftedJetProducerByMatchedObject.cc
Go to the documentation of this file.
2 
8 
9 template <typename T>
11  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
12 {
13  srcJets_ = consumes<JetCollection>(cfg.getParameter<edm::InputTag>("srcJets"));
14  srcUnshiftedObjects_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("srcUnshiftedObjects"));
15  srcShiftedObjects_ = consumes<edm::View<reco::Candidate> >(cfg.getParameter<edm::InputTag>("srcShiftedObjects"));
16 
17  dRmatch_Jet_ = cfg.getParameter<double>("dRmatch_Jet");
18  dRmatch_Object_ = cfg.exists("dRmatch_Object") ?
19  cfg.getParameter<double>("dRmatch_Object") : 0.1;
20 
22  dR2match_Object_ = dRmatch_Object_*dRmatch_Object_;
23 
24  produces<JetCollection>();
25 }
26 
27 template <typename T>
29 {
30 // nothing to be done yet...
31 }
32 
33 template <typename T>
35 {
36  edm::Handle<JetCollection> originalJets;
37  evt.getByToken(srcJets_, originalJets);
38 
39  edm::Handle<reco::CandidateView> unshiftedObjects;
40  evt.getByToken(srcUnshiftedObjects_, unshiftedObjects);
41 
42  edm::Handle<reco::CandidateView> shiftedObjects;
43  evt.getByToken(srcShiftedObjects_, shiftedObjects);
44 
45  objects_.clear();
46 
47  std::vector<bool> match(shiftedObjects->size(), false);
48  int prevMatch=-1;
49  int cnt = 0;
50 
51  for ( reco::CandidateView::const_iterator unshiftedObject = unshiftedObjects->begin();
52  unshiftedObject != unshiftedObjects->end(); ++unshiftedObject ) {
53  bool isMatched_Object = false;
54  double dR2bestMatch_Object = std::numeric_limits<double>::max();
55  prevMatch=-1;
56  cnt = 0;
57 
58  reco::Candidate::LorentzVector shiftedObjectP4_matched;
59  for ( reco::CandidateView::const_iterator shiftedObject = shiftedObjects->begin();
60  shiftedObject != shiftedObjects->end(); ++shiftedObject ) {
61  if( match[ cnt ] ) continue;
62 
63  double dR2 = deltaR2(unshiftedObject->p4(), shiftedObject->p4());
64  if ( dR2 < dR2match_Object_ && dR2 < dR2bestMatch_Object ) {
65  shiftedObjectP4_matched = shiftedObject->p4();
66  isMatched_Object = true;
67  dR2bestMatch_Object = dR2;
68 
69  prevMatch = cnt;
70  }
71  cnt++;
72  }
73  if ( isMatched_Object ) {
74  //Ambiguity removal
75  match[ prevMatch ] = true;
76  objects_.push_back(objectEntryType(shiftedObjectP4_matched, unshiftedObject->p4(), sqrt(dR2bestMatch_Object)));
77  }
78  }
79 
80 
81  match.assign(objects_.size(), false);
82 
83  std::auto_ptr<JetCollection> shiftedJets(new JetCollection);
84 
85  for ( typename JetCollection::const_iterator originalJet = originalJets->begin();
86  originalJet != originalJets->end(); ++originalJet ) {
87 
88  double shift = 0.;
89  bool applyShift = false;
90  double dR2bestMatch_Jet = std::numeric_limits<double>::max();
91  prevMatch=-1;
92  cnt = 0;
93 
94  for ( typename std::vector<objectEntryType>::const_iterator object = objects_.begin();
95  object != objects_.end(); ++object ) {
96  if ( !object->isValidMatch_ ) continue;
97  if( match[ cnt ] ) continue;
98 
99  double dR2 = deltaR2(originalJet->p4(), object->unshiftedObjectP4_);
100  if ( dR2 < dR2match_Jet_ && dR2 < dR2bestMatch_Jet ) {
101  shift = object->shift_;
102  applyShift = true;
103  dR2bestMatch_Jet = dR2;
104 
105  prevMatch = cnt;
106  }
107  cnt++;
108  }
109 
110  reco::Candidate::LorentzVector shiftedJetP4 = originalJet->p4();
111  if ( applyShift ) {
112  //Ambiguity removal
113  match[ prevMatch ] = true;
114 
115  shiftedJetP4 *= (1. + shift);
116  }
117 
118  T shiftedJet(*originalJet);
119  shiftedJet.setP4(shiftedJetP4);
120 
121  shiftedJets->push_back(shiftedJet);
122  }
123 
124  evt.put(shiftedJets);
125 }
126 
129 
132 
134 
137 
138 
139 
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:259
ShiftedJetProducerByMatchedObjectT< reco::CaloJet > ShiftedCaloJetProducerByMatchedObject
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
void produce(edm::Event &, const edm::EventSetup &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
T sqrt(T t)
Definition: SSEVec.h:48
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
edm::EDGetTokenT< edm::View< reco::Candidate > > srcShiftedObjects_
ShiftedJetProducerByMatchedObjectT< reco::PFJet > ShiftedPFJetProducerByMatchedObject
edm::EDGetTokenT< edm::View< reco::Candidate > > srcUnshiftedObjects_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
list object
Definition: dbtoconf.py:77
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85
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
moduleLabel_(iConfig.getParameter< string >("@module_label"))
long double T
ShiftedJetProducerByMatchedObjectT(const edm::ParameterSet &)