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  usePuppi_(iConfig.getParameter<bool>("usePuppi")) {
29 
30  if (residualsFromTxt_) {
31  residualsTxt_ = iConfig.getParameter<edm::FileInPath>("residualsTxt");
32  }
33 
34  std::vector<edm::ParameterSet> algos = iConfig.getParameter<std::vector<edm::ParameterSet> >("algos");
35  for (auto const& algoPset : algos) {
36  vAlgoGBRForestsAndConstants_.emplace_back(algoPset, runMvas_);
37  }
38 
39  if (!runMvas_) {
40  assert(algos.size() == 1);
41  }
42 }
43 
44 // ------------------------------------------------------------------------------------------
46 {
47  if ( globalCache->produceJetIds() ) {
48  produces<edm::ValueMap<StoredPileupJetIdentifier> > ("");
49  }
50  for (auto const& algoGBRForestsAndConstants : globalCache->vAlgoGBRForestsAndConstants()) {
51  std::string const& label = algoGBRForestsAndConstants.label();
52  algos_.emplace_back(label, std::make_unique<PileupJetIdAlgo>(&algoGBRForestsAndConstants));
53  if ( globalCache->runMvas() ) {
54  produces<edm::ValueMap<float> > (label+"Discriminant");
55  produces<edm::ValueMap<int> > (label+"Id");
56  }
57  }
58 
59  input_jet_token_ = consumes<edm::View<reco::Jet> >(iConfig.getParameter<edm::InputTag>("jets"));
60  input_vertex_token_ = consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertexes"));
61  input_vm_pujetid_token_ = consumes<edm::ValueMap<StoredPileupJetIdentifier> >(iConfig.getParameter<edm::InputTag>("jetids"));
62  input_rho_token_ = consumes<double>(iConfig.getParameter<edm::InputTag>("rho"));
63 }
64 
65 
66 
67 // ------------------------------------------------------------------------------------------
69 {
70 }
71 
72 
73 // ------------------------------------------------------------------------------------------
74 void
76 {
77  GBRForestsAndConstants const* gc = globalCache();
78 
79  using namespace edm;
80  using namespace std;
81  using namespace reco;
82 
83  // Input jets
84  Handle<View<Jet> > jetHandle;
85  iEvent.getByToken(input_jet_token_,jetHandle);
86  const View<Jet> & jets = *jetHandle;
87  // vertexes
88  Handle<VertexCollection> vertexHandle;
89  if( gc->produceJetIds() ) {
90  iEvent.getByToken(input_vertex_token_, vertexHandle);
91  }
92  const VertexCollection & vertexes = *(vertexHandle.product());
93  // input variables
95  if( ! gc->produceJetIds() ) {
96  iEvent.getByToken(input_vm_pujetid_token_, vmap);
97  }
98  // rho
100  double rho = 0.;
101 
102  // products
103  vector<StoredPileupJetIdentifier> ids;
104  map<string, vector<float> > mvas;
105  map<string, vector<int> > idflags;
106 
107  VertexCollection::const_iterator vtx;
108  if( gc->produceJetIds() ) {
109  // require basic quality cuts on the vertexes
110  vtx = vertexes.begin();
111  while( vtx != vertexes.end() && ( vtx->isFake() || vtx->ndof() < 4 ) ) {
112  ++vtx;
113  }
114  if( vtx == vertexes.end() ) { vtx = vertexes.begin(); }
115  }
116 
117  // Loop over input jets
118  bool ispat = true;
119  for ( unsigned int i=0; i<jets.size(); ++i ) {
120  // Pick the first algo to compute the input variables
121  auto algoi = algos_.begin();
122  PileupJetIdAlgo * ialgo = algoi->second.get();
123 
124  const Jet & jet = jets.at(i);
125  const pat::Jet * patjet = nullptr;
126  if(ispat) {
127  patjet=dynamic_cast<const pat::Jet *>(&jet);
128  ispat = patjet != nullptr;
129  }
130 
131  // Get jet energy correction
132  float jec = 0.;
133  if( gc->applyJec() ) {
134  // If haven't done it get rho from the event
135  if( rho == 0. ) {
136  iEvent.getByToken(input_rho_token_,rhoH);
137  rho = *rhoH;
138  }
139  // jet corrector
140  if( not jecCor_ ) {
141  initJetEnergyCorrector( iSetup, iEvent.isRealData() );
142  }
143  if( ispat ) {
144  jecCor_->setJetPt(patjet->correctedJet(0).pt());
145  } else {
146  jecCor_->setJetPt(jet.pt());
147  }
148  jecCor_->setJetEta(jet.eta());
149  jecCor_->setJetA(jet.jetArea());
150  jecCor_->setRho(rho);
151  jec = jecCor_->getCorrection();
152  }
153  // If it was requested AND the input is an uncorrected jet apply the JEC
154  bool applyJec = gc->applyJec() && ( ispat || !gc->inputIsCorrected() );
155  std::unique_ptr<reco::Jet> corrJet;
156 
157  if( applyJec ) {
158  float scale = jec;
159  if( ispat ) {
160  corrJet.reset(new pat::Jet(patjet->correctedJet(0))) ;
161  } else {
162  corrJet.reset(dynamic_cast<reco::Jet *>( jet.clone() ));
163  }
164  corrJet->scaleEnergy(scale);
165  }
166  const reco::Jet * theJet = ( applyJec ? corrJet.get() : &jet );
167 
168  PileupJetIdentifier puIdentifier;
169  if( gc->produceJetIds() ) {
170  // Compute the input variables
171  puIdentifier = ialgo->computeIdVariables(theJet, jec, &(*vtx), vertexes, rho,gc->usePuppi());
172  ids.push_back( puIdentifier );
173  } else {
174  // Or read it from the value map
175  puIdentifier = (*vmap)[jets.refAt(i)];
176  puIdentifier.jetPt(theJet->pt()); // make sure JEC is applied when computing the MVA
177  puIdentifier.jetEta(theJet->eta());
178  puIdentifier.jetPhi(theJet->phi());
179  ialgo->set(puIdentifier);
180  puIdentifier = ialgo->computeMva();
181  }
182 
183  if( gc->runMvas() ) {
184  // Compute the MVA and WP
185  mvas[algoi->first].push_back( puIdentifier.mva() );
186  idflags[algoi->first].push_back( puIdentifier.idFlag() );
187  for( ++algoi; algoi!=algos_.end(); ++algoi) {
188  ialgo = algoi->second.get();
189  ialgo->set(puIdentifier);
190  PileupJetIdentifier id = ialgo->computeMva();
191  mvas[algoi->first].push_back( id.mva() );
192  idflags[algoi->first].push_back( id.idFlag() );
193  }
194  }
195  }
196 
197  // Produce the output value maps
198  if( gc->runMvas() ) {
199  for( const auto& ialgo : algos_) {
200  // MVA
201  vector<float> & mva = mvas[ialgo.first];
202  auto mvaout = std::make_unique<ValueMap<float>>();
203  ValueMap<float>::Filler mvafiller(*mvaout);
204  mvafiller.insert(jetHandle,mva.begin(),mva.end());
205  mvafiller.fill();
206  iEvent.put(std::move(mvaout),ialgo.first+"Discriminant");
207 
208  // WP
209  vector<int> & idflag = idflags[ialgo.first];
210  auto idflagout = std::make_unique<ValueMap<int>>();
211  ValueMap<int>::Filler idflagfiller(*idflagout);
212  idflagfiller.insert(jetHandle,idflag.begin(),idflag.end());
213  idflagfiller.fill();
214  iEvent.put(std::move(idflagout),ialgo.first+"Id");
215  }
216  }
217  // input variables
218  if( gc->produceJetIds() ) {
219  assert( jetHandle->size() == ids.size() );
220  auto idsout = std::make_unique<ValueMap<StoredPileupJetIdentifier>>();
222  idsfiller.insert(jetHandle,ids.begin(),ids.end());
223  idsfiller.fill();
224  iEvent.put(std::move(idsout));
225  }
226 }
227 
228 
229 
230 // ------------------------------------------------------------------------------------------
231 void
233  //The following says we do not know what parameters are allowed so do no validation
234  // Please change this to state exactly what you do use, even if it is no parameters
236  desc.setUnknown();
237  descriptions.addDefault(desc);
238 }
239 
240 
241 // ------------------------------------------------------------------------------------------
242 void
244 {
245  GBRForestsAndConstants const* gc = globalCache();
246 
247  //jet energy correction levels to apply on raw jet
248  std::vector<std::string> jecLevels;
249  jecLevels.push_back("L1FastJet");
250  jecLevels.push_back("L2Relative");
251  jecLevels.push_back("L3Absolute");
252  if(isData && ! gc->residualsFromTxt() ) jecLevels.push_back("L2L3Residual");
253 
254  //check the corrector parameters needed according to the correction levels
256  iSetup.get<JetCorrectionsRecord>().get(gc->jec(),parameters);
257  for(std::vector<std::string>::const_iterator ll = jecLevels.begin(); ll != jecLevels.end(); ++ll)
258  {
259  const JetCorrectorParameters& ip = (*parameters)[*ll];
260  jetCorPars_.push_back(ip);
261  }
262  if( isData && gc->residualsFromTxt() ) {
263  jetCorPars_.push_back(JetCorrectorParameters(gc->residualsTxt().fullPath()));
264  }
265 
266  //instantiate the jet corrector
267  jecCor_ = std::make_unique<FactorizedJetCorrector>(jetCorPars_);
268 }
269 //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)
char const * label
PileupJetIdentifier computeIdVariables(const reco::Jet *jet, float jec, const reco::Vertex *, const reco::VertexCollection &, double rho, bool usePuppi)
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_
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:68
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:511
const float & jetEta() const