CMS 3D CMS Logo

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::global::EDProducer<> edm::global::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 QGTagger (const edm::ParameterSet &)
 
 ~QGTagger () override
 
- Public Member Functions inherited from edm::global::EDProducer<>
 EDProducer ()=default
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
bool wantsStreamLuminosityBlocks () const final
 
bool wantsStreamRuns () const final
 
- Public Member Functions inherited from edm::global::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 Descriptions method. More...
 
- Static Public Member Functions inherited from edm::global::EDProducerBase
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 > &, bool) const
 Calculation of axis2, mult and ptD. More...
 
bool isPackedCandidate (const reco::Jet *jet) const
 Function to tell us if we are using packedCandidates, only test for first candidate. More...
 
void produce (edm::StreamID, edm::Event &, const edm::EventSetup &) const override
 Produce qgLikelihood using {mult, ptD, -log(axis2)}. More...
 
template<typename T >
void putInEvent (const std::string &, const edm::Handle< edm::View< reco::Jet >> &, std::vector< T > *, edm::Event &) const
 Function to put product into event. More...
 

Private Attributes

edm::EDGetTokenT< reco::JetCorrectorjetCorrectorToken
 
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::VertexCollectionvertexToken
 

Additional Inherited Members

- Public Types inherited from edm::global::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
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 ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
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 16 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 28 of file QGTagger.cc.

References produceSyst, and qgLikelihood.

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

Member Function Documentation

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

Calculation of axis2, mult and ptD.

Definition at line 138 of file QGTagger.cc.

References a, b, EnergyCorrector::c, allConversions_cfi::d0, delta, reco::deltaPhi(), PVValHelper::dz, reco::LeafCandidate::eta(), reco::Jet::getJetConstituentsQuick(), edm::Ref< C, T, F >::isNonnull(), reco::LeafCandidate::phi(), funct::pow(), jets_cff::ptD, reco::TrackBase::qualityByName(), mathSSE::sqrt(), useQC, extraflags_cff::vtx, and mps_merge::weight.

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

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

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

Referenced by ~QGTagger().

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

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

Definition at line 127 of file QGTagger.cc.

References Exception, and reco::Jet::getJetConstituentsQuick().

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

127  {
128  for(auto candidate : jet->getJetConstituentsQuick()){
129  if(typeid(pat::PackedCandidate)==typeid(*candidate)) return true;
130  else if(typeid(reco::PFCandidate)==typeid(*candidate)) return false;
131  else throw cms::Exception("WrongJetCollection", "Jet constituents are not particle flow candidates");
132  }
133  return false;
134 }
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
virtual std::vector< const reco::Candidate * > getJetConstituentsQuick() const
quick list of constituents
void QGTagger::produce ( edm::StreamID  ,
edm::Event iEvent,
const edm::EventSetup iSetup 
) const
overrideprivate

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

Definition at line 53 of file QGTagger.cc.

References calcVariables(), QGLikelihoodCalculator::computeQGLikelihood(), reco::JetCorrector::correction(), objects.autophobj::float, edm::EventSetup::get(), edm::eventsetup::EventSetupRecordImplementation< T >::get(), edm::Event::getByToken(), isPackedCandidate(), metsig::jet, jetCorrectorToken, fwrapper::jets, jetsLabel, jetsToken, cmsBatch::log, produceSyst, EnergyCorrector::pt, jets_cff::ptD, putInEvent(), qgLikelihood, rho, rhoToken, QGLikelihoodCalculator::systematicSmearing(), systLabel, useJetCorr, particleFlowSuperClusterECAL_cfi::vertexCollection, and vertexToken.

Referenced by ~QGTagger().

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

Function to put product into event.

Definition at line 116 of file QGTagger.cc.

References objects.autophobj::filler, fwrapper::jets, eostools::move(), MillePedeFileConverter_cfg::out, and edm::Event::put().

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

116  {
117  auto out = std::make_unique<edm::ValueMap<T>>();
119  filler.insert(jets, product->begin(), product->end());
120  filler.fill();
121  iEvent.put(std::move(out), name);
122  delete product;
123 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
def move(src, dest)
Definition: eostools.py:511

Member Data Documentation

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

Definition at line 29 of file QGTagger.h.

Referenced by produce().

std::string QGTagger::jetsLabel
private

Definition at line 32 of file QGTagger.h.

Referenced by produce().

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

Definition at line 28 of file QGTagger.h.

Referenced by produce().

const bool QGTagger::produceSyst
private

Definition at line 33 of file QGTagger.h.

Referenced by produce(), and QGTagger().

QGLikelihoodCalculator* QGTagger::qgLikelihood
private

Definition at line 34 of file QGTagger.h.

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

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

Definition at line 31 of file QGTagger.h.

Referenced by produce().

std::string QGTagger::systLabel
private

Definition at line 32 of file QGTagger.h.

Referenced by produce().

const bool QGTagger::useJetCorr
private

Definition at line 33 of file QGTagger.h.

Referenced by produce().

const bool QGTagger::useQC
private

Definition at line 33 of file QGTagger.h.

Referenced by calcVariables().

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

Definition at line 30 of file QGTagger.h.

Referenced by produce().