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 
20 
22  runMvas_(iConfig.getParameter<bool>("runMvas")),
23  produceJetIds_(iConfig.getParameter<bool>("produceJetIds")),
24  inputIsCorrected_(iConfig.getParameter<bool>("inputIsCorrected")),
25  applyJec_(iConfig.getParameter<bool>("applyJec")),
26  jec_(iConfig.getParameter<std::string>("jec")),
27  residualsFromTxt_(iConfig.getParameter<bool>("residualsFromTxt")) {
28 
29  if (residualsFromTxt_) {
30  residualsTxt_ = iConfig.getParameter<edm::FileInPath>("residualsTxt");
31  }
32 
33  std::vector<edm::ParameterSet> algos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
34  for (auto const& algoPset : algos) {
35  vAlgoGBRForestsAndConstants_.emplace_back(algoPset, runMvas_);
36  }
37 
38  if (!runMvas_) {
39  assert(algos.size() == 1);
40  }
41 }
42 
43 // ------------------------------------------------------------------------------------------
45 {
46  if ( globalCache->produceJetIds() ) {
47  produces<edm::ValueMap<StoredPileupJetIdentifier> > ("");
48  }
49  for (auto const& algoGBRForestsAndConstants : globalCache->vAlgoGBRForestsAndConstants()) {
50  std::string const& label = algoGBRForestsAndConstants.label();
51  algos_.emplace_back(label, std::make_unique<PileupJetIdAlgo>(&algoGBRForestsAndConstants));
52  if ( globalCache->runMvas() ) {
53  produces<edm::ValueMap<float> > (label+"Discriminant");
54  produces<edm::ValueMap<int> > (label+"Id");
55  }
56  }
57 
58  input_jet_token_ = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
59  input_vertex_token_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexes"));
60  input_vm_pujetid_token_ = consumes<edm::ValueMap<StoredPileupJetIdentifier> >(iConfig.getParameter<edm::InputTag>("jetids"));
61  input_rho_token_ = consumes<double>(iConfig.getParameter<edm::InputTag>("rho"));
62 }
63 
64 
65 
66 // ------------------------------------------------------------------------------------------
68 {
69 }
70 
71 
72 // ------------------------------------------------------------------------------------------
73 void
75 {
76  GBRForestsAndConstants const* gc = globalCache();
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( gc->produceJetIds() ) {
89  iEvent.getByToken(input_vertex_token_, vertexHandle);
90  }
91  const VertexCollection & vertexes = *(vertexHandle.product());
92  // input variables
94  if( ! gc->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( gc->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  bool ispat = true;
118  for ( unsigned int i=0; i<jets.size(); ++i ) {
119  // Pick the first algo to compute the input variables
120  auto algoi = algos_.begin();
121  PileupJetIdAlgo * ialgo = algoi->second.get();
122 
123  const Jet & jet = jets.at(i);
124  const pat::Jet * patjet = nullptr;
125  if(ispat) {
126  patjet=dynamic_cast<const pat::Jet *>(&jet);
127  ispat = patjet != nullptr;
128  }
129 
130  // Get jet energy correction
131  float jec = 0.;
132  if( gc->applyJec() ) {
133  // If haven't done it get rho from the event
134  if( rho == 0. ) {
135  iEvent.getByToken(input_rho_token_,rhoH);
136  rho = *rhoH;
137  }
138  // jet corrector
139  if( not jecCor_ ) {
140  initJetEnergyCorrector( iSetup, iEvent.isRealData() );
141  }
142  if( ispat ) {
143  jecCor_->setJetPt(patjet->correctedJet(0).pt());
144  } else {
145  jecCor_->setJetPt(jet.pt());
146  }
147  jecCor_->setJetEta(jet.eta());
148  jecCor_->setJetA(jet.jetArea());
149  jecCor_->setRho(rho);
150  jec = jecCor_->getCorrection();
151  }
152  // If it was requested AND the input is an uncorrected jet apply the JEC
153  bool applyJec = gc->applyJec() && ( ispat || !gc->inputIsCorrected() );
154  std::unique_ptr<reco::Jet> corrJet;
155 
156  if( applyJec ) {
157  float scale = jec;
158  if( ispat ) {
159  corrJet.reset(new pat::Jet(patjet->correctedJet(0))) ;
160  } else {
161  corrJet.reset(dynamic_cast<reco::Jet *>( jet.clone() ));
162  }
163  corrJet->scaleEnergy(scale);
164  }
165  const reco::Jet * theJet = ( applyJec ? corrJet.get() : &jet );
166 
167  PileupJetIdentifier puIdentifier;
168  if( gc->produceJetIds() ) {
169  // Compute the input variables
170  puIdentifier = ialgo->computeIdVariables(theJet, jec, &(*vtx), vertexes, rho);
171  ids.push_back( puIdentifier );
172  } else {
173  // Or read it from the value map
174  puIdentifier = (*vmap)[jets.refAt(i)];
175  puIdentifier.jetPt(theJet->pt()); // make sure JEC is applied when computing the MVA
176  puIdentifier.jetEta(theJet->eta());
177  puIdentifier.jetPhi(theJet->phi());
178  ialgo->set(puIdentifier);
179  puIdentifier = ialgo->computeMva();
180  }
181 
182  if( gc->runMvas() ) {
183  // Compute the MVA and WP
184  mvas[algoi->first].push_back( puIdentifier.mva() );
185  idflags[algoi->first].push_back( puIdentifier.idFlag() );
186  for( ++algoi; algoi!=algos_.end(); ++algoi) {
187  ialgo = algoi->second.get();
188  ialgo->set(puIdentifier);
189  PileupJetIdentifier id = ialgo->computeMva();
190  mvas[algoi->first].push_back( id.mva() );
191  idflags[algoi->first].push_back( id.idFlag() );
192  }
193  }
194  }
195 
196  // Produce the output value maps
197  if( gc->runMvas() ) {
198  for( const auto& ialgo : algos_) {
199  // MVA
200  vector<float> & mva = mvas[ialgo.first];
201  auto mvaout = std::make_unique<ValueMap<float>>();
202  ValueMap<float>::Filler mvafiller(*mvaout);
203  mvafiller.insert(jetHandle,mva.begin(),mva.end());
204  mvafiller.fill();
205  iEvent.put(std::move(mvaout),ialgo.first+"Discriminant");
206 
207  // WP
208  vector<int> & idflag = idflags[ialgo.first];
209  auto idflagout = std::make_unique<ValueMap<int>>();
210  ValueMap<int>::Filler idflagfiller(*idflagout);
211  idflagfiller.insert(jetHandle,idflag.begin(),idflag.end());
212  idflagfiller.fill();
213  iEvent.put(std::move(idflagout),ialgo.first+"Id");
214  }
215  }
216  // input variables
217  if( gc->produceJetIds() ) {
218  assert( jetHandle->size() == ids.size() );
219  auto idsout = std::make_unique<ValueMap<StoredPileupJetIdentifier>>();
221  idsfiller.insert(jetHandle,ids.begin(),ids.end());
222  idsfiller.fill();
223  iEvent.put(std::move(idsout));
224  }
225 }
226 
227 
228 
229 // ------------------------------------------------------------------------------------------
230 void
232  //The following says we do not know what parameters are allowed so do no validation
233  // Please change this to state exactly what you do use, even if it is no parameters
235  desc.setUnknown();
236  descriptions.addDefault(desc);
237 }
238 
239 
240 // ------------------------------------------------------------------------------------------
241 void
243 {
244  GBRForestsAndConstants const* gc = globalCache();
245 
246  //jet energy correction levels to apply on raw jet
247  std::vector<std::string> jecLevels;
248  jecLevels.push_back("L1FastJet");
249  jecLevels.push_back("L2Relative");
250  jecLevels.push_back("L3Absolute");
251  if(isData && ! gc->residualsFromTxt() ) jecLevels.push_back("L2L3Residual");
252 
253  //check the corrector parameters needed according to the correction levels
255  iSetup.get<JetCorrectionsRecord>().get(gc->jec(),parameters);
256  for(std::vector<std::string>::const_iterator ll = jecLevels.begin(); ll != jecLevels.end(); ++ll)
257  {
258  const JetCorrectorParameters& ip = (*parameters)[*ll];
259  jetCorPars_.push_back(ip);
260  }
261  if( isData && gc->residualsFromTxt() ) {
262  jetCorPars_.push_back(JetCorrectorParameters(gc->residualsTxt().fullPath()));
263  }
264 
265  //instantiate the jet corrector
266  jecCor_ = std::make_unique<FactorizedJetCorrector>(jetCorPars_);
267 }
268 //define this as a plug-in
std::vector< PileupJetIdAlgo::AlgoGBRForestsAndConstants > const & vAlgoGBRForestsAndConstants() const
T getParameter(std::string const &) const
void set(const PileupJetIdentifier &)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
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
size_type size() const
double pt() const final
transverse momentum
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 &, GBRForestsAndConstants const *)
RefToBase< value_type > refAt(size_type i) const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
Definition: Jet.py:1
std::string const & jec() const
vector< PseudoJet > jets
edm::FileInPath const & residualsTxt() const
void initJetEnergyCorrector(const edm::EventSetup &iSetup, bool isData)
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< PileupJetIdAlgo::AlgoGBRForestsAndConstants > vAlgoGBRForestsAndConstants_
const float & jetPhi() const
edm::FileInPath residualsTxt_
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho)
T const * product() const
Definition: Handle.h:81
PileupJetIdentifier computeMva()
Analysis-level calorimeter jet class.
Definition: Jet.h:80
CompositePtrCandidate * clone() const override
returns a clone of the candidate
fixed size matrix
HLT enums.
GBRForestsAndConstants(edm::ParameterSet const &)
virtual float jetArea() const
get jet area
Definition: Jet.h:105
T get() const
Definition: EventSetup.h:63
const_reference at(size_type pos) const
std::string fullPath() const
Definition: FileInPath.cc:197
Jet correctedJet(const std::string &level, const std::string &flavor="none", const std::string &set="") const
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:510
const float & jetEta() const