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 //
17 //
18 
19 #include <memory>
21 
22 
23 // ------------------------------------------------------------------------------------------
25 {
26  runMvas_ = iConfig.getParameter<bool>("runMvas");
27  produceJetIds_ = iConfig.getParameter<bool>("produceJetIds");
28  jets_ = iConfig.getParameter<edm::InputTag>("jets");
29  vertexes_ = iConfig.getParameter<edm::InputTag>("vertexes");
30  jetids_ = iConfig.getParameter<edm::InputTag>("jetids");
31  inputIsCorrected_ = iConfig.getParameter<bool>("inputIsCorrected");
32  applyJec_ = iConfig.getParameter<bool>("applyJec");
33  jec_ = iConfig.getParameter<std::string>("jec");
34  rho_ = iConfig.getParameter<edm::InputTag>("rho");
35  residualsFromTxt_ = iConfig.getParameter<bool>("residualsFromTxt");
36  if(residualsFromTxt_) residualsTxt_ = iConfig.getParameter<edm::FileInPath>("residualsTxt");
37  std::vector<edm::ParameterSet> algos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
38 
39  jecCor_ = 0;
40 
41  if( ! runMvas_ ) assert( algos.size() == 1 );
42 
43  if( produceJetIds_ ) {
44  produces<edm::ValueMap<StoredPileupJetIdentifier> > ("");
45  }
46  for(std::vector<edm::ParameterSet>::iterator it=algos.begin(); it!=algos.end(); ++it) {
47  std::string label = it->getParameter<std::string>("label");
48  algos_.push_back( std::make_pair(label,new PileupJetIdAlgo(*it, runMvas_)) );
49  if( runMvas_ ) {
50  produces<edm::ValueMap<float> > (label+"Discriminant");
51  produces<edm::ValueMap<int> > (label+"Id");
52  }
53  }
54 
55  input_jet_token_ = consumes<edm::View<reco::Jet> >(jets_);
56  input_vertex_token_ = consumes<reco::VertexCollection>(vertexes_);
57  input_vm_pujetid_token_ = consumes<edm::ValueMap<StoredPileupJetIdentifier> >(jetids_);
58  input_rho_token_ = consumes<double>(rho_);
59 
60 }
61 
62 
63 
64 // ------------------------------------------------------------------------------------------
66 {
67  for(std::vector<std::pair<std::string,PileupJetIdAlgo *> >::iterator
68  ialgo = algos_.begin(); ialgo!=algos_.end(); ++ialgo) {
69  delete ialgo->second;
70  }
71 }
72 
73 
74 // ------------------------------------------------------------------------------------------
75 void
77 {
78  using namespace edm;
79  using namespace std;
80  using namespace reco;
81 
82  // Input jets
83  Handle<View<Jet> > jetHandle;
84  iEvent.getByToken(input_jet_token_,jetHandle);
85  const View<Jet> & jets = *jetHandle;
86  // vertexes
87  Handle<VertexCollection> vertexHandle;
88  if( produceJetIds_ ) {
89  iEvent.getByToken(input_vertex_token_, vertexHandle);
90  }
91  const VertexCollection & vertexes = *(vertexHandle.product());
92  // input variables
94  if( ! produceJetIds_ ) {
95  iEvent.getByToken(input_vm_pujetid_token_, vmap);
96  }
97  // rho
99  double rho = 0.;
100 
101  // products
102  vector<StoredPileupJetIdentifier> ids;
103  map<string, vector<float> > mvas;
104  map<string, vector<int> > idflags;
105 
106  VertexCollection::const_iterator vtx;
107  if( produceJetIds_ ) {
108  // require basic quality cuts on the vertexes
109  vtx = vertexes.begin();
110  while( vtx != vertexes.end() && ( vtx->isFake() || vtx->ndof() < 4 ) ) {
111  ++vtx;
112  }
113  if( vtx == vertexes.end() ) { vtx = vertexes.begin(); }
114  }
115 
116  // Loop over input jets
117  for ( unsigned int i=0; i<jets.size(); ++i ) {
118  // Pick the first algo to compute the input variables
119  vector<pair<string,PileupJetIdAlgo *> >::iterator algoi = algos_.begin();
120  PileupJetIdAlgo * ialgo = algoi->second;
121 
122  const Jet & jet = jets.at(i);
123  //const pat::Jet * patjet = dynamic_cast<const pat::Jet *>(&jet);
124  //bool ispat = patjet != 0;
125 
126  // Get jet energy correction
127  float jec = 0.;
128  if( applyJec_ ) {
129  // If haven't done it get rho from the event
130  if( rho == 0. ) {
131  iEvent.getByToken(input_rho_token_,rhoH);
132  rho = *rhoH;
133  }
134  // jet corrector
135  if( jecCor_ == 0 ) {
136  initJetEnergyCorrector( iSetup, iEvent.isRealData() );
137  }
138  //if( ispat ) {
139  // jecCor_->setJetPt(patjet->correctedJet(0).pt());
140  //} else {
141  jecCor_->setJetPt(jet.pt());
142  //}
143  jecCor_->setJetEta(jet.eta());
144  jecCor_->setJetA(jet.jetArea());
145  jecCor_->setRho(rho);
146  jec = jecCor_->getCorrection();
147  }
148  // If it was requested AND the input is an uncorrected jet apply the JEC
149  bool applyJec = applyJec_ && !inputIsCorrected_; //( ! ispat && ! inputIsCorrected_ );
150  reco::Jet * corrJet = 0;
151 
152  if( applyJec ) {
153  float scale = jec;
154  //if( ispat ) {
155  // corrJet = new pat::Jet(patjet->correctedJet(0)) ;
156  //} else {
157  corrJet = dynamic_cast<reco::Jet *>( jet.clone() );
158  //}
159  corrJet->scaleEnergy(scale);
160  }
161  const reco::Jet * theJet = ( applyJec ? corrJet : &jet );
162 
163  PileupJetIdentifier puIdentifier;
164  if( produceJetIds_ ) {
165  // Compute the input variables
166  puIdentifier = ialgo->computeIdVariables(theJet, jec, &(*vtx), vertexes, rho);
167  ids.push_back( puIdentifier );
168  } else {
169  // Or read it from the value map
170  puIdentifier = (*vmap)[jets.refAt(i)];
171  puIdentifier.jetPt(theJet->pt()); // make sure JEC is applied when computing the MVA
172  puIdentifier.jetEta(theJet->eta());
173  puIdentifier.jetPhi(theJet->phi());
174  ialgo->set(puIdentifier);
175  puIdentifier = ialgo->computeMva();
176  }
177 
178  if( runMvas_ ) {
179  // Compute the MVA and WP
180  mvas[algoi->first].push_back( puIdentifier.mva() );
181  idflags[algoi->first].push_back( puIdentifier.idFlag() );
182  for( ++algoi; algoi!=algos_.end(); ++algoi) {
183  ialgo = algoi->second;
184  ialgo->set(puIdentifier);
185  PileupJetIdentifier id = ialgo->computeMva();
186  mvas[algoi->first].push_back( id.mva() );
187  idflags[algoi->first].push_back( id.idFlag() );
188  }
189  }
190 
191  // cleanup
192  if( corrJet ) { delete corrJet; }
193  }
194 
195  // Produce the output value maps
196  if( runMvas_ ) {
197  for(vector<pair<string,PileupJetIdAlgo *> >::iterator ialgo = algos_.begin(); ialgo!=algos_.end(); ++ialgo) {
198  // MVA
199  vector<float> & mva = mvas[ialgo->first];
200  auto_ptr<ValueMap<float> > mvaout(new ValueMap<float>());
201  ValueMap<float>::Filler mvafiller(*mvaout);
202  mvafiller.insert(jetHandle,mva.begin(),mva.end());
203  mvafiller.fill();
204  iEvent.put(mvaout,ialgo->first+"Discriminant");
205 
206  // WP
207  vector<int> & idflag = idflags[ialgo->first];
208  auto_ptr<ValueMap<int> > idflagout(new ValueMap<int>());
209  ValueMap<int>::Filler idflagfiller(*idflagout);
210  idflagfiller.insert(jetHandle,idflag.begin(),idflag.end());
211  idflagfiller.fill();
212  iEvent.put(idflagout,ialgo->first+"Id");
213  }
214  }
215  // input variables
216  if( produceJetIds_ ) {
217  assert( jetHandle->size() == ids.size() );
218  auto_ptr<ValueMap<StoredPileupJetIdentifier> > idsout(new ValueMap<StoredPileupJetIdentifier>());
220  idsfiller.insert(jetHandle,ids.begin(),ids.end());
221  idsfiller.fill();
222  iEvent.put(idsout);
223  }
224 }
225 
226 
227 
228 // ------------------------------------------------------------------------------------------
229 void
231  //The following says we do not know what parameters are allowed so do no validation
232  // Please change this to state exactly what you do use, even if it is no parameters
234  desc.setUnknown();
235  descriptions.addDefault(desc);
236 }
237 
238 
239 // ------------------------------------------------------------------------------------------
240 void
242 {
243  //jet energy correction levels to apply on raw jet
244  std::vector<std::string> jecLevels;
245  jecLevels.push_back("L1FastJet");
246  jecLevels.push_back("L2Relative");
247  jecLevels.push_back("L3Absolute");
248  if(isData && ! residualsFromTxt_ ) jecLevels.push_back("L2L3Residual");
249 
250  //check the corrector parameters needed according to the correction levels
252  iSetup.get<JetCorrectionsRecord>().get(jec_,parameters);
253  for(std::vector<std::string>::const_iterator ll = jecLevels.begin(); ll != jecLevels.end(); ++ll)
254  {
255  const JetCorrectorParameters& ip = (*parameters)[*ll];
256  jetCorPars_.push_back(ip);
257  }
258  if( isData && residualsFromTxt_ ) {
260  }
261 
262  //instantiate the jet corrector
264 }
265 //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_
virtual void scaleEnergy(double fScale)
scale energy of the jet
edm::EDGetTokenT< reco::VertexCollection > input_vertex_token_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
const int & idFlag() const
std::vector< std::pair< std::string, PileupJetIdAlgo * > > algos_
const float & mva() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Base class for all types of Jets.
Definition: Jet.h:20
assert(m_qm.get())
edm::EDGetTokenT< edm::View< reco::Jet > > input_jet_token_
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
bool isRealData() const
Definition: EventBase.h:64
const float & jetPt() const
PileupJetIdProducer(const edm::ParameterSet &)
virtual double eta() const
momentum pseudorapidity
virtual double pt() const
transverse momentum
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
virtual CompositePtrCandidate * clone() const
returns a clone of the candidate
FactorizedJetCorrector * jecCor_
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
vector< PseudoJet > jets
edm::EDGetTokenT< double > input_rho_token_
edm::FileInPath residualsTxt_
void initJetEnergyCorrector(const edm::EventSetup &iSetup, bool isData)
virtual void produce(edm::Event &, const edm::EventSetup &) override
const float & jetPhi() const
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho)
const T & get() const
Definition: EventSetup.h:56
PileupJetIdentifier computeMva()
virtual float jetArea() const
get jet area
Definition: Jet.h:105
edm::EDGetTokenT< edm::ValueMap< StoredPileupJetIdentifier > > input_vm_pujetid_token_
std::string fullPath() const
Definition: FileInPath.cc:165
virtual double phi() const
momentum azimuthal angle
const float & jetEta() const