CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
20 
24 
30 
31 template < class T >
33 
34 public:
35 
37 
39  srcToken_( consumes< typename edm::View<T> >( params.getParameter<edm::InputTag>("src") ) ),
40  matchedToken_( consumes< typename edm::View<T> >( 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_, StringObjectFunction<T>( value_, lazyParser_ ) ) );
51  produces< JetValueMap >();
52  }
53 
54  if( valueLabels_.size()>0 || values_.size()>0 )
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], StringObjectFunction<T>( 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  virtual void beginJob() override {}
75  virtual void endJob() override {}
76 
77  virtual void produce(edm::Event& iEvent, const edm::EventSetup& iSetup) override {
78 
80  iEvent.getByToken( srcToken_, h_jets1 );
82  iEvent.getByToken( matchedToken_, h_jets2 );
83 
84  std::vector<float> values( h_jets1->size(), -99999 );
85  std::map<std::string, std::vector<float> > valuesMap;
86  if( multiValue_ )
87  {
88  for( size_t i=0; i<valueLabels_.size(); ++i)
89  valuesMap.insert( std::make_pair( valueLabels_[i], std::vector<float>( h_jets1->size(), -99999 ) ) );
90  }
91  std::vector<bool> jets1_locks( h_jets1->size(), false );
92 
93  for ( typename edm::View<T>::const_iterator ibegin = h_jets2->begin(),
94  iend = h_jets2->end(), ijet = ibegin;
95  ijet != iend; ++ijet )
96  {
97  float matched_dR2 = 1e9;
98  int matched_index = -1;
99 
100  for ( typename edm::View<T>::const_iterator jbegin = h_jets1->begin(),
101  jend = h_jets1->end(), jjet = jbegin;
102  jjet != jend; ++jjet )
103  {
104  int index=jjet - jbegin;
105 
106  if( jets1_locks.at(index) ) continue; // skip jets that have already been matched
107 
108  float temp_dR2 = reco::deltaR2(ijet->eta(),ijet->phi(),jjet->eta(),jjet->phi());
109  if ( temp_dR2 < matched_dR2 )
110  {
111  matched_dR2 = temp_dR2;
112  matched_index = index;
113  }
114  }// end loop over src jets
115 
116  if( matched_index>=0 )
117  {
118  if ( matched_dR2 > distMax_*distMax_ )
119  edm::LogWarning("MatchedJetsFarApart") << "Matched jets separated by dR greater than distMax=" << distMax_;
120  else
121  {
122  jets1_locks.at(matched_index) = true;
123  if( value_!="" )
124  values.at(matched_index) = (evaluationMap_.at(value_))(*ijet);
125  if( multiValue_ )
126  {
127  for( size_t i=0; i<valueLabels_.size(); ++i)
128  valuesMap.at(valueLabels_[i]).at(matched_index) = (evaluationMap_.at(valueLabels_[i]))(*ijet);
129  }
130  }
131  }
132  }// end loop over matched jets
133 
134  if( value_!="" )
135  {
136  std::auto_ptr< JetValueMap > jetValueMap ( new JetValueMap() );
137 
138  JetValueMap::Filler filler(*jetValueMap);
139  filler.insert(h_jets1, values.begin(), values.end());
140  filler.fill();
141 
142  // put in Event
143  iEvent.put(jetValueMap);
144  }
145  if( multiValue_ )
146  {
147  for( size_t i=0; i<valueLabels_.size(); ++i)
148  {
149  std::auto_ptr< JetValueMap > jetValueMap ( new JetValueMap() );
150 
151  JetValueMap::Filler filler(*jetValueMap);
152  filler.insert(h_jets1, valuesMap.at(valueLabels_[i]).begin(), valuesMap.at(valueLabels_[i]).end());
153  filler.fill();
154 
155  // put in Event
156  iEvent.put(jetValueMap, valueLabels_[i]);
157  }
158  }
159  }
160 
163  const double distMax_;
165  const std::vector<std::string> values_;
166  const std::vector<std::string> valueLabels_;
167  const bool lazyParser_;
169  std::map<std::string, StringObjectFunction<T> > evaluationMap_;
170 };
171 
173 
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
const std::vector< std::string > values_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::ValueMap< float > JetValueMap
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
const edm::EDGetTokenT< typename edm::View< T > > srcToken_
const std::vector< std::string > valueLabels_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EDGetTokenT< typename edm::View< T > > matchedToken_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
std::map< std::string, StringObjectFunction< T > > evaluationMap_
virtual void beginJob() override
JetDeltaRValueMapProducer< reco::Jet > RecoJetDeltaRValueMapProducer
volatile std::atomic< bool > shutdown_flag false
long double T
JetDeltaRValueMapProducer(edm::ParameterSet const &params)