CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QGTagger.cc
Go to the documentation of this file.
1 #include <memory>
2 
10 
18 
23 
30  jetsToken( consumes<edm::View<reco::Jet>>( iConfig.getParameter<edm::InputTag>("srcJets"))),
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 }
52 
53 
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 }
110 
112 template <typename T> void QGTagger::putInEvent(std::string name, const edm::Handle<edm::View<reco::Jet>>& jets, std::vector<T>* product, edm::Event& iEvent){
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 }
120 
121 
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 }
132 
133 
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 }
211 
212 
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 }
225 
226 
227 //define this as a plug-in
dbl * delta
Definition: mlp_gen.cc:36
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:252
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Base class for all types of Jets.
Definition: Jet.h:20
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
edm::EDGetTokenT< reco::VertexCollection > vertexToken
Definition: QGTagger.h:28
const bool useJetCorr
Definition: QGTagger.h:31
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
virtual void produce(edm::Event &, const edm::EventSetup &)
Produce qgLikelihood using {mult, ptD, -log(axis2)}.
Definition: QGTagger.cc:55
QGLikelihoodCalculator * qgLikelihood
Definition: QGTagger.h:33
std::string systLabel
Definition: QGTagger.h:30
tuple vertexCollection
QGTagger(const edm::ParameterSet &)
Definition: QGTagger.cc:29
virtual double eta() const
momentum pseudorapidity
int iEvent
Definition: GenABIO.cc:230
const bool produceSyst
Definition: QGTagger.h:31
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
T sqrt(T t)
Definition: SSEVec.h:48
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
bool isPackedCandidate(const reco::Candidate *candidate)
Function to tell us if we are using packedCandidates, only test for first candidate.
Definition: QGTagger.cc:123
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Descriptions method.
Definition: QGTagger.cc:214
std::tuple< int, float, float > calcVariables(const reco::Jet *, edm::Handle< reco::VertexCollection > &)
Calculation of axis2, mult and ptD.
Definition: QGTagger.cc:135
tuple out
Definition: dbtoconf.py:99
edm::EDGetTokenT< edm::View< reco::Jet > > jetsToken
Definition: QGTagger.h:26
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
std::string jetsLabel
Definition: QGTagger.h:30
const T & get() const
Definition: EventSetup.h:56
double b
Definition: hdecay.h:120
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< double > rhoToken
Definition: QGTagger.h:29
bool weAreUsingPackedCandidates
Definition: QGTagger.h:32
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
double a
Definition: hdecay.h:121
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
bool weStillNeedToCheckJetCandidates
Definition: QGTagger.h:32
edm::EDGetTokenT< reco::JetCorrector > jetCorrectorToken
Definition: QGTagger.h:27
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
tuple log
Definition: cmsBatch.py:341
float systematicSmearing(edm::ESHandle< QGLikelihoodSystematicsObject > &QGLParamsColl, float pt, float eta, float rho, float qgValue, int qgIndex)