CMS 3D CMS Logo

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