CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes
QGTagger Class Reference

#include <QGTagger.h>

Inheritance diagram for QGTagger:
edm::EDProducer edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 QGTagger (const edm::ParameterSet &)
 
 ~QGTagger ()
 
- Public Member Functions inherited from edm::EDProducer
 EDProducer ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDProducer ()
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 Descriptions method. More...
 
- Static Public Member Functions inherited from edm::EDProducer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 

Private Member Functions

std::tuple< int, float, float > calcVariables (const reco::Jet *, edm::Handle< reco::VertexCollection > &)
 Calculation of axis2, mult and ptD. More...
 
bool isPackedCandidate (const reco::Candidate *candidate)
 Function to tell us if we are using packedCandidates, only test for first candidate. More...
 
virtual void produce (edm::Event &, const edm::EventSetup &)
 Produce qgLikelihood using {mult, ptD, -log(axis2)}. More...
 
template<typename T >
void putInEvent (std::string, const edm::Handle< edm::View< reco::Jet >> &, std::vector< T > *, edm::Event &)
 Function to put product into event. More...
 

Private Attributes

edm::EDGetTokenT
< reco::JetCorrector
jetCorrectorToken
 
std::string jetsLabel
 
edm::EDGetTokenT< edm::View
< reco::Jet > > 
jetsToken
 
const bool produceSyst
 
QGLikelihoodCalculatorqgLikelihood
 
edm::EDGetTokenT< double > rhoToken
 
std::string systLabel
 
const bool useJetCorr
 
const bool useQC
 
edm::EDGetTokenT
< reco::VertexCollection
vertexToken
 
bool weAreUsingPackedCandidates
 
bool weStillNeedToCheckJetCandidates
 

Additional Inherited Members

- Public Types inherited from edm::EDProducer
typedef EDProducer ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 14 of file QGTagger.h.

Constructor & Destructor Documentation

QGTagger::QGTagger ( const edm::ParameterSet iConfig)
explicit

EDProducer class to produced the qgLikelihood values and related variables If the input jets are uncorrected, the jecService should be provided, so jet are corrected on the fly before the algorithm is applied Authors: andre.nosp@m.a.ca.nosp@m.rlo.m.nosp@m.arin.nosp@m.i@cer.nosp@m.n.ch, tom.c.nosp@m.orne.nosp@m.lis@c.nosp@m.ern..nosp@m.ch, cms-q.nosp@m.g-wo.nosp@m.rking.nosp@m.grou.nosp@m.p@cer.nosp@m.n.ch

Definition at line 29 of file QGTagger.cc.

References produceSyst, qgLikelihood, and weStillNeedToCheckJetCandidates.

29  :
31  jetCorrectorToken( consumes<reco::JetCorrector>( iConfig.getParameter<edm::InputTag>("jec"))),
32  vertexToken( consumes<reco::VertexCollection>( iConfig.getParameter<edm::InputTag>("srcVertexCollection"))),
33  rhoToken( consumes<double>( iConfig.getParameter<edm::InputTag>("srcRho"))),
34  jetsLabel( iConfig.getParameter<std::string>("jetsLabel")),
35  systLabel( iConfig.getParameter<std::string>("systematicsLabel")),
36  useQC( iConfig.getParameter<bool>("useQualityCuts")),
37  useJetCorr( !iConfig.getParameter<edm::InputTag>("jec").label().empty()),
38  produceSyst( systLabel != "")
39 {
40  produces<edm::ValueMap<float>>("qgLikelihood");
41  produces<edm::ValueMap<float>>("axis2");
42  produces<edm::ValueMap<int>>("mult");
43  produces<edm::ValueMap<float>>("ptD");
44  if(produceSyst){
45  produces<edm::ValueMap<float>>("qgLikelihoodSmearedQuark");
46  produces<edm::ValueMap<float>>("qgLikelihoodSmearedGluon");
47  produces<edm::ValueMap<float>>("qgLikelihoodSmearedAll");
48  }
51 }
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::VertexCollection > vertexToken
Definition: QGTagger.h:28
const bool useJetCorr
Definition: QGTagger.h:31
QGLikelihoodCalculator * qgLikelihood
Definition: QGTagger.h:33
std::string systLabel
Definition: QGTagger.h:30
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const bool produceSyst
Definition: QGTagger.h:31
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken
Definition: QGTagger.h:26
std::string jetsLabel
Definition: QGTagger.h:30
edm::EDGetTokenT< double > rhoToken
Definition: QGTagger.h:29
std::string const & label() const
Definition: InputTag.h:43
bool weStillNeedToCheckJetCandidates
Definition: QGTagger.h:32
edm::EDGetTokenT< reco::JetCorrector > jetCorrectorToken
Definition: QGTagger.h:27
const bool useQC
Definition: QGTagger.h:31
QGTagger::~QGTagger ( )
inline

Definition at line 17 of file QGTagger.h.

References qgLikelihood.

17 { delete qgLikelihood;};
QGLikelihoodCalculator * qgLikelihood
Definition: QGTagger.h:33

Member Function Documentation

std::tuple< int, float, float > QGTagger::calcVariables ( const reco::Jet jet,
edm::Handle< reco::VertexCollection > &  vC 
)
private

Calculation of axis2, mult and ptD.

Definition at line 135 of file QGTagger.cc.

References a, b, EnergyCorrector::c, delta, reco::deltaPhi(), reco::LeafCandidate::eta(), reco::Jet::getJetConstituentsQuick(), edm::Ref< C, T, F >::isNonnull(), isPackedCandidate(), VarParsing::mult, reco::LeafCandidate::phi(), funct::pow(), reco::TrackBase::qualityByName(), mathSSE::sqrt(), useQC, and puppiForMET_cff::weight.

Referenced by produce().

135  {
136  float sum_weight = 0., sum_deta = 0., sum_dphi = 0., sum_deta2 = 0., sum_dphi2 = 0., sum_detadphi = 0., sum_pt = 0.;
137  int mult = 0;
138 
139  //Loop over the jet constituents
140  for(auto daughter : jet->getJetConstituentsQuick()){
141  if(isPackedCandidate(daughter)){ //packed candidate situation
142  auto part = static_cast<const pat::PackedCandidate*>(daughter);
143 
144  if(part->charge()){
145  if(!(part->fromPV() > 1 && part->trackHighPurity())) continue;
146  if(useQC){
147  if((part->dz()*part->dz())/(part->dzError()*part->dzError()) > 25.) continue;
148  if((part->dxy()*part->dxy())/(part->dxyError()*part->dxyError()) < 25.) ++mult;
149  } else ++mult;
150  } else {
151  if(part->pt() < 1.0) continue;
152  ++mult;
153  }
154  } else {
155  auto part = static_cast<const reco::PFCandidate*>(daughter);
156 
157  reco::TrackRef itrk = part->trackRef();
158  if(itrk.isNonnull()){ //Track exists --> charged particle
159  auto vtxLead = vC->begin();
160  auto vtxClose = vC->begin(); //Search for closest vertex to track
161  for(auto vtx = vC->begin(); vtx != vC->end(); ++vtx){
162  if(fabs(itrk->dz(vtx->position())) < fabs(itrk->dz(vtxClose->position()))) vtxClose = vtx;
163  }
164  if(!(vtxClose == vtxLead && itrk->quality(reco::TrackBase::qualityByName("highPurity")))) continue;
165 
166  if(useQC){ //If useQC, require dz and d0 cuts
167  float dz = itrk->dz(vtxClose->position());
168  float d0 = itrk->dxy(vtxClose->position());
169  float dz_sigma_square = pow(itrk->dzError(),2) + pow(vtxClose->zError(),2);
170  float d0_sigma_square = pow(itrk->d0Error(),2) + pow(vtxClose->xError(),2) + pow(vtxClose->yError(),2);
171  if(dz*dz/dz_sigma_square > 25.) continue;
172  if(d0*d0/d0_sigma_square < 25.) ++mult;
173  } else ++mult;
174  } else { //No track --> neutral particle
175  if(part->pt() < 1.0) continue; //Only use neutrals with pt > 1 GeV
176  ++mult;
177  }
178  }
179 
180  float deta = daughter->eta() - jet->eta();
181  float dphi = reco::deltaPhi(daughter->phi(), jet->phi());
182  float partPt = daughter->pt();
183  float weight = partPt*partPt;
184 
185  sum_weight += weight;
186  sum_pt += partPt;
187  sum_deta += deta*weight;
188  sum_dphi += dphi*weight;
189  sum_deta2 += deta*deta*weight;
190  sum_detadphi += deta*dphi*weight;
191  sum_dphi2 += dphi*dphi*weight;
192  }
193 
194  //Calculate axis2 and ptD
195  float a = 0., b = 0., c = 0.;
196  float ave_deta = 0., ave_dphi = 0., ave_deta2 = 0., ave_dphi2 = 0.;
197  if(sum_weight > 0){
198  ave_deta = sum_deta/sum_weight;
199  ave_dphi = sum_dphi/sum_weight;
200  ave_deta2 = sum_deta2/sum_weight;
201  ave_dphi2 = sum_dphi2/sum_weight;
202  a = ave_deta2 - ave_deta*ave_deta;
203  b = ave_dphi2 - ave_dphi*ave_dphi;
204  c = -(sum_detadphi/sum_weight - ave_deta*ave_dphi);
205  }
206  float delta = sqrt(fabs((a-b)*(a-b)+4*c*c));
207  float axis2 = (a+b-delta > 0 ? sqrt(0.5*(a+b-delta)) : 0);
208  float ptD = (sum_weight > 0 ? sqrt(sum_weight)/sum_pt : 0);
209  return std::make_tuple(mult, ptD, axis2);
210 }
dbl * delta
Definition: mlp_gen.cc:36
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
virtual double eta() const
momentum pseudorapidity
T sqrt(T t)
Definition: SSEVec.h:48
bool isPackedCandidate(const reco::Candidate *candidate)
Function to tell us if we are using packedCandidates, only test for first candidate.
Definition: QGTagger.cc:123
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:123
part
Definition: HCALResponse.h:20
double b
Definition: hdecay.h:120
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
double a
Definition: hdecay.h:121
virtual std::vector< const reco::Candidate * > getJetConstituentsQuick() const
quick list of constituents
virtual double phi() const
momentum azimuthal angle
const bool useQC
Definition: QGTagger.h:31
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
void QGTagger::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Descriptions method.

Definition at line 214 of file QGTagger.cc.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), HLT_25ns14e33_v1_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

214  {
216  desc.add<edm::InputTag>("srcJets");
217  desc.add<edm::InputTag>("srcRho");
218  desc.add<std::string>("jetsLabel");
219  desc.add<std::string>("systematicsLabel", "");
220  desc.add<bool>("useQualityCuts");
221  desc.add<edm::InputTag>("jec", edm::InputTag())->setComment("Jet correction service: only applied when non-empty");
222  desc.add<edm::InputTag>("srcVertexCollection")->setComment("Ignored for miniAOD, possible to keep empty");
223  descriptions.add("QGTagger", desc);
224 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool QGTagger::isPackedCandidate ( const reco::Candidate candidate)
private

Function to tell us if we are using packedCandidates, only test for first candidate.

Definition at line 123 of file QGTagger.cc.

References Exception, weAreUsingPackedCandidates, and weStillNeedToCheckJetCandidates.

Referenced by calcVariables().

123  {
125  if(typeid(pat::PackedCandidate)==typeid(*candidate)) weAreUsingPackedCandidates = true;
126  else if(typeid(reco::PFCandidate)==typeid(*candidate)) weAreUsingPackedCandidates = false;
127  else throw cms::Exception("WrongJetCollection", "Jet constituents are not particle flow candidates");
129  }
131 }
bool weAreUsingPackedCandidates
Definition: QGTagger.h:32
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
bool weStillNeedToCheckJetCandidates
Definition: QGTagger.h:32
void QGTagger::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
)
privatevirtual

Produce qgLikelihood using {mult, ptD, -log(axis2)}.

Implements edm::EDProducer.

Definition at line 55 of file QGTagger.cc.

References calcVariables(), QGLikelihoodCalculator::computeQGLikelihood(), edm::EventSetup::get(), edm::eventsetup::EventSetupRecord::get(), edm::Event::getByToken(), metsig::jet, jetCorrectorToken, fwrapper::jets, jetsLabel, jetsToken, cmsBatch::log, VarParsing::mult, produceSyst, EnergyCorrector::pt, putInEvent(), qgLikelihood, rho, rhoToken, QGLikelihoodCalculator::systematicSmearing(), systLabel, useJetCorr, GoodVertex_cfg::vertexCollection, and vertexToken.

55  {
56  std::vector<float>* qgProduct = new std::vector<float>;
57  std::vector<float>* axis2Product = new std::vector<float>;
58  std::vector<int>* multProduct = new std::vector<int>;
59  std::vector<float>* ptDProduct = new std::vector<float>;
60  std::vector<float>* smearedQuarkProduct = new std::vector<float>;
61  std::vector<float>* smearedGluonProduct = new std::vector<float>;
62  std::vector<float>* smearedAllProduct = new std::vector<float>;
63 
68 
70  QGLikelihoodRcd const & rcdhandle = iSetup.get<QGLikelihoodRcd>();
71  rcdhandle.get(jetsLabel, QGLParamsColl);
72 
74  if(produceSyst){
75  QGLikelihoodSystematicsRcd const & systrcdhandle = iSetup.get<QGLikelihoodSystematicsRcd>();
76  systrcdhandle.get(systLabel, QGLSystColl);
77  }
78 
79  for(auto jet = jets->begin(); jet != jets->end(); ++jet){
80  float pt = (useJetCorr ? jet->pt()*jetCorr->correction(*jet) : jet->pt());
81 
82  float ptD, axis2; int mult;
83  std::tie(mult, ptD, axis2) = calcVariables(&*jet, vertexCollection);
84 
85  float qgValue;
86  if(mult > 2) qgValue = qgLikelihood->computeQGLikelihood(QGLParamsColl, pt, jet->eta(), *rho, {(float) mult, ptD, -std::log(axis2)});
87  else qgValue = -1;
88 
89  qgProduct->push_back(qgValue);
90  if(produceSyst){
91  smearedQuarkProduct->push_back(qgLikelihood->systematicSmearing(QGLSystColl, pt, jet->eta(), *rho, qgValue, 0));
92  smearedGluonProduct->push_back(qgLikelihood->systematicSmearing(QGLSystColl, pt, jet->eta(), *rho, qgValue, 1));
93  smearedAllProduct->push_back(qgLikelihood->systematicSmearing( QGLSystColl, pt, jet->eta(), *rho, qgValue, 2));
94  }
95  axis2Product->push_back(axis2);
96  multProduct->push_back(mult);
97  ptDProduct->push_back(ptD);
98  }
99 
100  putInEvent("qgLikelihood", jets, qgProduct, iEvent);
101  putInEvent("axis2", jets, axis2Product, iEvent);
102  putInEvent("mult", jets, multProduct, iEvent);
103  putInEvent("ptD", jets, ptDProduct, iEvent);
104  if(produceSyst){
105  putInEvent("qgLikelihoodSmearedQuark", jets, smearedQuarkProduct, iEvent);
106  putInEvent("qgLikelihoodSmearedGluon", jets, smearedGluonProduct, iEvent);
107  putInEvent("qgLikelihoodSmearedAll", jets, smearedAllProduct, iEvent);
108  }
109 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
edm::EDGetTokenT< reco::VertexCollection > vertexToken
Definition: QGTagger.h:28
const bool useJetCorr
Definition: QGTagger.h:31
QGLikelihoodCalculator * qgLikelihood
Definition: QGTagger.h:33
std::string systLabel
Definition: QGTagger.h:30
tuple vertexCollection
const bool produceSyst
Definition: QGTagger.h:31
float computeQGLikelihood(edm::ESHandle< QGLikelihoodObject > &QGLParamsColl, float pt, float eta, float rho, std::vector< float > vars)
Compute likelihood for a jet using the QGLikelihoodObject information and a set of variables...
vector< PseudoJet > jets
void get(HolderT &iHolder) const
std::tuple< int, float, float > calcVariables(const reco::Jet *, edm::Handle< reco::VertexCollection > &)
Calculation of axis2, mult and ptD.
Definition: QGTagger.cc:135
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken
Definition: QGTagger.h:26
std::string jetsLabel
Definition: QGTagger.h:30
const T & get() const
Definition: EventSetup.h:56
edm::EDGetTokenT< double > rhoToken
Definition: QGTagger.h:29
void putInEvent(std::string, const edm::Handle< edm::View< reco::Jet >> &, std::vector< T > *, edm::Event &)
Function to put product into event.
Definition: QGTagger.cc:112
edm::EDGetTokenT< reco::JetCorrector > jetCorrectorToken
Definition: QGTagger.h:27
tuple log
Definition: cmsBatch.py:341
float systematicSmearing(edm::ESHandle< QGLikelihoodSystematicsObject > &QGLParamsColl, float pt, float eta, float rho, float qgValue, int qgIndex)
template<typename T >
void QGTagger::putInEvent ( std::string  name,
const edm::Handle< edm::View< reco::Jet >> &  jets,
std::vector< T > *  product,
edm::Event iEvent 
)
private

Function to put product into event.

Definition at line 112 of file QGTagger.cc.

References edm::helper::Filler< Map >::fill(), edm::helper::Filler< Map >::insert(), fwrapper::jets, dbtoconf::out, and edm::Event::put().

Referenced by produce().

112  {
113  std::auto_ptr<edm::ValueMap<T>> out(new edm::ValueMap<T>());
114  typename edm::ValueMap<T>::Filler filler(*out);
115  filler.insert(jets, product->begin(), product->end());
116  filler.fill();
117  iEvent.put(out, name);
118  delete product;
119 }
const_iterator begin() const
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
tuple out
Definition: dbtoconf.py:99

Member Data Documentation

edm::EDGetTokenT<reco::JetCorrector> QGTagger::jetCorrectorToken
private

Definition at line 27 of file QGTagger.h.

Referenced by produce().

std::string QGTagger::jetsLabel
private

Definition at line 30 of file QGTagger.h.

Referenced by produce().

edm::EDGetTokenT<edm::View<reco::Jet> > QGTagger::jetsToken
private

Definition at line 26 of file QGTagger.h.

Referenced by produce().

const bool QGTagger::produceSyst
private

Definition at line 31 of file QGTagger.h.

Referenced by produce(), and QGTagger().

QGLikelihoodCalculator* QGTagger::qgLikelihood
private

Definition at line 33 of file QGTagger.h.

Referenced by produce(), QGTagger(), and ~QGTagger().

edm::EDGetTokenT<double> QGTagger::rhoToken
private

Definition at line 29 of file QGTagger.h.

Referenced by produce().

std::string QGTagger::systLabel
private

Definition at line 30 of file QGTagger.h.

Referenced by produce().

const bool QGTagger::useJetCorr
private

Definition at line 31 of file QGTagger.h.

Referenced by produce().

const bool QGTagger::useQC
private

Definition at line 31 of file QGTagger.h.

Referenced by calcVariables().

edm::EDGetTokenT<reco::VertexCollection> QGTagger::vertexToken
private

Definition at line 28 of file QGTagger.h.

Referenced by produce().

bool QGTagger::weAreUsingPackedCandidates
private

Definition at line 32 of file QGTagger.h.

Referenced by isPackedCandidate().

bool QGTagger::weStillNeedToCheckJetCandidates
private

Definition at line 32 of file QGTagger.h.

Referenced by isPackedCandidate(), and QGTagger().