CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PileupJetIdProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PileupJetIdProducer
4 // Class: PileupJetIdProducer
5 //
13 //
14 // Original Author: Pasquale Musella,40 2-A12,+41227671706,
15 // Created: Wed Apr 18 15:48:47 CEST 2012
16 // $Id: PileupJetIdProducer.cc,v 1.3 2013/02/27 20:57:37 eulisse Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
28 
34 
36 
40 
44 
46 
47 // ------------------------------------------------------------------------------------------
49 public:
50  explicit PileupJetIdProducer(const edm::ParameterSet&);
52 
53  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
54 
55 private:
56  virtual void produce(edm::Event&, const edm::EventSetup&);
57 
58 
59  void initJetEnergyCorrector(const edm::EventSetup &iSetup, bool isData);
60 
64  std::vector<std::pair<std::string, PileupJetIdAlgo *> > algos_;
65 
69  std::vector<JetCorrectorParameters> jetCorPars_;
70 };
71 
72 // ------------------------------------------------------------------------------------------
74 {
75  runMvas_ = iConfig.getParameter<bool>("runMvas");
76  produceJetIds_ = iConfig.getParameter<bool>("produceJetIds");
77  jets_ = iConfig.getParameter<edm::InputTag>("jets");
78  vertexes_ = iConfig.getParameter<edm::InputTag>("vertexes");
79  jetids_ = iConfig.getParameter<edm::InputTag>("jetids");
80  inputIsCorrected_ = iConfig.getParameter<bool>("inputIsCorrected");
81  applyJec_ = iConfig.getParameter<bool>("applyJec");
82  jec_ = iConfig.getParameter<std::string>("jec");
83  rho_ = iConfig.getParameter<edm::InputTag>("rho");
84  residualsFromTxt_ = iConfig.getParameter<bool>("residualsFromTxt");
85  residualsTxt_ = iConfig.getParameter<edm::FileInPath>("residualsTxt");
86  std::vector<edm::ParameterSet> algos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
87 
88  jecCor_ = 0;
89 
90  if( ! runMvas_ ) assert( algos.size() == 1 );
91 
92  if( produceJetIds_ ) {
93  produces<edm::ValueMap<StoredPileupJetIdentifier> > ("");
94  }
95  for(std::vector<edm::ParameterSet>::iterator it=algos.begin(); it!=algos.end(); ++it) {
96  std::string label = it->getParameter<std::string>("label");
97  algos_.push_back( std::make_pair(label,new PileupJetIdAlgo(*it)) );
98  if( runMvas_ ) {
99  produces<edm::ValueMap<float> > (label+"Discriminant");
100  produces<edm::ValueMap<int> > (label+"Id");
101  }
102  }
103 }
104 
105 
106 
107 // ------------------------------------------------------------------------------------------
109 {
110 }
111 
112 
113 // ------------------------------------------------------------------------------------------
114 void
116 {
117  using namespace edm;
118  using namespace std;
119  using namespace reco;
120 
121  // Input jets
122  Handle<View<Jet> > jetHandle;
123  iEvent.getByLabel(jets_,jetHandle);
124  const View<Jet> & jets = *jetHandle;
125  // vertexes
126  Handle<VertexCollection> vertexHandle;
127  if( produceJetIds_ ) {
128  iEvent.getByLabel(vertexes_, vertexHandle);
129  }
130  const VertexCollection & vertexes = *(vertexHandle.product());
131  // input variables
133  if( ! produceJetIds_ ) {
134  iEvent.getByLabel(jetids_, vmap);
135  }
136  // rho
138  double rho = 0.;
139 
140  // products
141  vector<StoredPileupJetIdentifier> ids;
142  map<string, vector<float> > mvas;
143  map<string, vector<int> > idflags;
144 
145  VertexCollection::const_iterator vtx;
146  if( produceJetIds_ ) {
147  // require basic quality cuts on the vertexes
148  vtx = vertexes.begin();
149  while( vtx != vertexes.end() && ( vtx->isFake() || vtx->ndof() < 4 ) ) {
150  ++vtx;
151  }
152  if( vtx == vertexes.end() ) { vtx = vertexes.begin(); }
153  }
154 
155  // Loop over input jets
156  for ( unsigned int i=0; i<jets.size(); ++i ) {
157  // Pick the first algo to compute the input variables
158  vector<pair<string,PileupJetIdAlgo *> >::iterator algoi = algos_.begin();
159  PileupJetIdAlgo * ialgo = algoi->second;
160 
161  const Jet & jet = jets.at(i);
162  //const pat::Jet * patjet = dynamic_cast<const pat::Jet *>(&jet);
163  //bool ispat = patjet != 0;
164 
165  // Get jet energy correction
166  float jec = 0.;
167  if( applyJec_ ) {
168  // If haven't done it get rho from the event
169  if( rho == 0. ) {
170  iEvent.getByLabel(rho_,rhoH);
171  rho = *rhoH;
172  }
173  // jet corrector
174  if( jecCor_ == 0 ) {
175  initJetEnergyCorrector( iSetup, iEvent.isRealData() );
176  }
177  //if( ispat ) {
178  // jecCor_->setJetPt(patjet->correctedJet(0).pt());
179  //} else {
180  jecCor_->setJetPt(jet.pt());
181  //}
182  jecCor_->setJetEta(jet.eta());
183  jecCor_->setJetA(jet.jetArea());
184  jecCor_->setRho(rho);
185  jec = jecCor_->getCorrection();
186  }
187 
188  // If it was requested or the input is an uncorrected jet apply the JEC
189  bool applyJec = applyJec_ || !inputIsCorrected_; //( ! ispat && ! inputIsCorrected_ );
190  reco::Jet * corrJet = 0;
191  if( applyJec ) {
192  float scale = jec;
193  //if( ispat ) {
194  // corrJet = new pat::Jet(patjet->correctedJet(0)) ;
195  //} else {
196  corrJet = dynamic_cast<reco::Jet *>( jet.clone() );
197  //}
198  corrJet->scaleEnergy(scale);
199  }
200  const reco::Jet * theJet = ( applyJec ? corrJet : &jet );
201 
202  PileupJetIdentifier puIdentifier;
203  if( produceJetIds_ ) {
204  // Compute the input variables
205  puIdentifier = ialgo->computeIdVariables(theJet, jec, &(*vtx), vertexes, runMvas_);
206  ids.push_back( puIdentifier );
207  } else {
208  // Or read it from the value map
209  puIdentifier = (*vmap)[jets.refAt(i)];
210  puIdentifier.jetPt(theJet->pt()); // make sure JEC is applied when computing the MVA
211  puIdentifier.jetEta(theJet->eta());
212  puIdentifier.jetPhi(theJet->phi());
213  ialgo->set(puIdentifier);
214  puIdentifier = ialgo->computeMva();
215  }
216 
217  if( runMvas_ ) {
218  // Compute the MVA and WP
219  mvas[algoi->first].push_back( puIdentifier.mva() );
220  idflags[algoi->first].push_back( puIdentifier.idFlag() );
221  for( ++algoi; algoi!=algos_.end(); ++algoi) {
222  ialgo = algoi->second;
223  ialgo->set(puIdentifier);
224  PileupJetIdentifier id = ialgo->computeMva();
225  mvas[algoi->first].push_back( id.mva() );
226  idflags[algoi->first].push_back( id.idFlag() );
227  }
228  }
229 
230  // cleanup
231  if( corrJet ) { delete corrJet; }
232  }
233 
234  // Produce the output value maps
235  if( runMvas_ ) {
236  for(vector<pair<string,PileupJetIdAlgo *> >::iterator ialgo = algos_.begin(); ialgo!=algos_.end(); ++ialgo) {
237  // MVA
238  vector<float> & mva = mvas[ialgo->first];
239  auto_ptr<ValueMap<float> > mvaout(new ValueMap<float>());
240  ValueMap<float>::Filler mvafiller(*mvaout);
241  mvafiller.insert(jetHandle,mva.begin(),mva.end());
242  mvafiller.fill();
243  iEvent.put(mvaout,ialgo->first+"Discriminant");
244 
245  // WP
246  vector<int> & idflag = idflags[ialgo->first];
247  auto_ptr<ValueMap<int> > idflagout(new ValueMap<int>());
248  ValueMap<int>::Filler idflagfiller(*idflagout);
249  idflagfiller.insert(jetHandle,idflag.begin(),idflag.end());
250  idflagfiller.fill();
251  iEvent.put(idflagout,ialgo->first+"Id");
252  }
253  }
254  // input variables
255  if( produceJetIds_ ) {
256  assert( jetHandle->size() == ids.size() );
257  auto_ptr<ValueMap<StoredPileupJetIdentifier> > idsout(new ValueMap<StoredPileupJetIdentifier>());
259  idsfiller.insert(jetHandle,ids.begin(),ids.end());
260  idsfiller.fill();
261  iEvent.put(idsout);
262  }
263 }
264 
265 
266 
267 // ------------------------------------------------------------------------------------------
268 void
270  //The following says we do not know what parameters are allowed so do no validation
271  // Please change this to state exactly what you do use, even if it is no parameters
273  desc.setUnknown();
274  descriptions.addDefault(desc);
275 }
276 
277 
278 // ------------------------------------------------------------------------------------------
279 void
281 {
282  //jet energy correction levels to apply on raw jet
283  std::vector<std::string> jecLevels;
284  jecLevels.push_back("L1FastJet");
285  jecLevels.push_back("L2Relative");
286  jecLevels.push_back("L3Absolute");
287  if(isData && ! residualsFromTxt_ ) jecLevels.push_back("L2L3Residual");
288 
289  //check the corrector parameters needed according to the correction levels
291  iSetup.get<JetCorrectionsRecord>().get(jec_,parameters);
292  for(std::vector<std::string>::const_iterator ll = jecLevels.begin(); ll != jecLevels.end(); ++ll)
293  {
294  const JetCorrectorParameters& ip = (*parameters)[*ll];
295  jetCorPars_.push_back(ip);
296  }
297  if( isData && residualsFromTxt_ ) {
299  }
300 
301  //instantiate the jet corrector
303 }
304 //define this as a plug-in
T getParameter(std::string const &) const
void set(const PileupJetIdentifier &)
int i
Definition: DBlmapReader.cc:9
dictionary parameters
Definition: Parameters.py:2
std::vector< JetCorrectorParameters > jetCorPars_
std::vector< std::pair< std::string, PileupJetIdAlgo * > > algos_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Base class for all types of Jets.
Definition: Jet.h:21
Definition: DDAxes.h:10
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
virtual void scaleEnergy(double fScale)
scale energy of the jet
Definition: Jet.cc:445
bool isRealData() const
Definition: EventBase.h:60
PileupJetIdProducer(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual CompositePtrCandidate * clone() const
returns a clone of the candidate
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
FactorizedJetCorrector * jecCor_
int iEvent
Definition: GenABIO.cc:243
void addDefault(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
vector< PseudoJet > jets
edm::FileInPath residualsTxt_
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, bool calculateMva=false)
void initJetEnergyCorrector(const edm::EventSetup &iSetup, bool isData)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
const T & get() const
Definition: EventSetup.h:55
PileupJetIdentifier computeMva()
virtual float jetArea() const
get jet area
Definition: Jet.h:106
std::string fullPath() const
Definition: FileInPath.cc:171
virtual float pt() const GCC11_FINAL
transverse momentum
virtual void produce(edm::Event &, const edm::EventSetup &)