CMS 3D CMS Logo

JetDeltaRValueMapProducer.cc
Go to the documentation of this file.
1 /* \class JetDeltaRValueMapProducer
2  *
3  * Associates jets using delta-R matching, and writes out
4  * a valuemap of single float variables based on a StringObjectFunction.
5  * This is used for miniAOD.
6  *
7  * \author: Sal Rappoccio
8  *
9  *
10  * for more details about the cut syntax, see the documentation
11  * page below:
12  *
13  * https://twiki.cern.ch/twiki/bin/view/CMS/SWGuidePhysicsCutParser
14  *
15  *
16  */
17 
18 
24 
30 
31 template < class T, class C=T >
33 
34 public:
35 
37 
39  srcToken_( consumes< typename edm::View<T> >( params.getParameter<edm::InputTag>("src") ) ),
40  matchedToken_( consumes< typename edm::View<C> >( params.getParameter<edm::InputTag>( "matched" ) ) ),
41  distMax_( params.getParameter<double>( "distMax" ) ),
42  value_( params.existsAs<std::string>("value") ? params.getParameter<std::string>("value") : "" ),
43  values_( params.existsAs<std::vector<std::string> >("values") ? params.getParameter<std::vector<std::string> >("values") : std::vector<std::string>() ),
44  valueLabels_( params.existsAs<std::vector<std::string> >("valueLabels") ? params.getParameter<std::vector<std::string> >("valueLabels") : std::vector<std::string>() ),
45  lazyParser_( params.existsAs<bool>("lazyParser") ? params.getParameter<bool>("lazyParser") : false ),
47  {
48  if( value_!="" )
49  {
50  evaluationMap_.insert( std::make_pair( value_, std::unique_ptr<StringObjectFunction<C> >( new StringObjectFunction<C>( value_, lazyParser_ ) ) ) );
51  produces< JetValueMap >();
52  }
53 
54  if( !valueLabels_.empty() || !values_.empty() )
55  {
56  if( valueLabels_.size()==values_.size() )
57  {
58  multiValue_ = true;
59  for( size_t i=0; i<valueLabels_.size(); ++i)
60  {
61  evaluationMap_.insert( std::make_pair( valueLabels_[i], std::unique_ptr<StringObjectFunction<C> >( new StringObjectFunction<C>( values_[i], lazyParser_ ) ) ) );
62  produces< JetValueMap >(valueLabels_[i]);
63  }
64  }
65  else
66  edm::LogWarning("ValueLabelMismatch") << "The number of value labels does not match the number of values. Values will not be evaluated.";
67  }
68  }
69 
71 
72 private:
73 
74  void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override {
75 
77  iEvent.getByToken( srcToken_, h_jets1 );
79  iEvent.getByToken( matchedToken_, h_jets2 );
80 
81  std::vector<float> values( h_jets1->size(), -99999 );
82  std::map<std::string, std::vector<float> > valuesMap;
83  if( multiValue_ )
84  {
85  for( size_t i=0; i<valueLabels_.size(); ++i)
86  valuesMap.insert( std::make_pair( valueLabels_[i], std::vector<float>( h_jets1->size(), -99999 ) ) );
87  }
88  std::vector<bool> jets1_locks( h_jets1->size(), false );
89 
90  for ( typename edm::View<C>::const_iterator ibegin = h_jets2->begin(),
91  iend = h_jets2->end(), ijet = ibegin;
92  ijet != iend; ++ijet )
93  {
94  float matched_dR2 = 1e9;
95  int matched_index = -1;
96 
97  for ( typename edm::View<T>::const_iterator jbegin = h_jets1->begin(),
98  jend = h_jets1->end(), jjet = jbegin;
99  jjet != jend; ++jjet )
100  {
101  int index=jjet - jbegin;
102 
103  if( jets1_locks.at(index) ) continue; // skip jets that have already been matched
104 
105  float temp_dR2 = reco::deltaR2(ijet->eta(),ijet->phi(),jjet->eta(),jjet->phi());
106  if ( temp_dR2 < matched_dR2 )
107  {
108  matched_dR2 = temp_dR2;
109  matched_index = index;
110  }
111  }// end loop over src jets
112 
113  if( matched_index>=0 )
114  {
115  if ( matched_dR2 > distMax_*distMax_ )
116  edm::LogInfo("MatchedJetsFarApart") << "Matched jets separated by dR greater than distMax=" << distMax_;
117  else
118  {
119  jets1_locks.at(matched_index) = true;
120  if( value_!="" )
121  values.at(matched_index) = (*(evaluationMap_.at(value_)))(*ijet);
122  if( multiValue_ )
123  {
124  for( size_t i=0; i<valueLabels_.size(); ++i)
125  valuesMap.at(valueLabels_[i]).at(matched_index) = (*(evaluationMap_.at(valueLabels_[i])))(*ijet);
126  }
127  }
128  }
129  }// end loop over matched jets
130 
131  if( value_!="" )
132  {
133  std::unique_ptr< JetValueMap > jetValueMap ( new JetValueMap() );
134 
135  JetValueMap::Filler filler(*jetValueMap);
136  filler.insert(h_jets1, values.begin(), values.end());
137  filler.fill();
138 
139  // put in Event
140  iEvent.put(std::move(jetValueMap));
141  }
142  if( multiValue_ )
143  {
144  for( size_t i=0; i<valueLabels_.size(); ++i)
145  {
146  std::unique_ptr< JetValueMap > jetValueMap ( new JetValueMap() );
147 
148  JetValueMap::Filler filler(*jetValueMap);
149  filler.insert(h_jets1, valuesMap.at(valueLabels_[i]).begin(), valuesMap.at(valueLabels_[i]).end());
150  filler.fill();
151 
152  // put in Event
153  iEvent.put(std::move(jetValueMap), valueLabels_[i]);
154  }
155  }
156  }
157 
160  const double distMax_;
162  const std::vector<std::string> values_;
163  const std::vector<std::string> valueLabels_;
164  const bool lazyParser_;
166  std::map<std::string, std::unique_ptr<const StringObjectFunction<C> > > evaluationMap_;
167 };
168 
172 
edm::ValueMap< float > JetValueMap
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
JetDeltaRValueMapProducer(edm::ParameterSet const &params)
const std::vector< std::string > valueLabels_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
JetDeltaRValueMapProducer< reco::Jet, pat::Jet > RecoJetToPatJetDeltaRValueMapProducer
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
std::map< std::string, std::unique_ptr< const StringObjectFunction< C > > > evaluationMap_
JetDeltaRValueMapProducer< pat::Jet > PatJetDeltaRValueMapProducer
const edm::EDGetTokenT< typename edm::View< C > > matchedToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:230
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
const edm::EDGetTokenT< typename edm::View< T > > srcToken_
T1 deltaR2(T1 eta1, T2 phi1, T3 eta2, T4 phi2)
Definition: deltaR.h:36
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
JetDeltaRValueMapProducer< reco::Jet > RecoJetDeltaRValueMapProducer
const std::vector< std::string > values_
long double T
def move(src, dest)
Definition: eostools.py:510