Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055 #include "FWCore/Framework/interface/EDProducer.h"
00056 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00057 #include "FWCore/ParameterSet/interface/ParameterSetfwd.h"
00058 #include "FWCore/Utilities/interface/InputTag.h"
00059
00060 #include "FWCore/Framework/interface/Event.h"
00061 #include "FWCore/Framework/interface/EventSetup.h"
00062 #include "FWCore/Framework/interface/MakerMacros.h"
00063 #include "FWCore/Framework/interface/ESHandle.h"
00064 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00065
00066 #include "DataFormats/JetReco/interface/Jet.h"
00067 #include "DataFormats/JetReco/interface/JetCollection.h"
00068 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00069 #include "SimDataFormats/JetMatching/interface/JetFlavour.h"
00070 #include "SimDataFormats/JetMatching/interface/JetFlavourMatching.h"
00071 #include "SimDataFormats/JetMatching/interface/MatchedPartons.h"
00072 #include "SimDataFormats/JetMatching/interface/JetMatchedPartons.h"
00073
00074 #include "DataFormats/Common/interface/View.h"
00075 #include "DataFormats/Common/interface/Ref.h"
00076 #include "DataFormats/Common/interface/getRef.h"
00077
00078
00079 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00080 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00081
00082 #include <memory>
00083 #include <string>
00084 #include <iostream>
00085 #include <vector>
00086 #include <Math/VectorUtil.h>
00087 #include <TMath.h>
00088
00089 using namespace std;
00090 using namespace reco;
00091 using namespace edm;
00092 using namespace ROOT::Math::VectorUtil;
00093
00094 class JetPartonMatcher : public edm::EDProducer
00095 {
00096 public:
00097 JetPartonMatcher( const edm::ParameterSet & );
00098 ~JetPartonMatcher();
00099
00100 private:
00101 virtual void produce(edm::Event&, const edm::EventSetup& );
00102
00103 int fillAlgoritDefinition( const Jet& );
00104 int fillPhysicsDefinition( const Jet& );
00105 int theHeaviest;
00106 int theNearest2;
00107 int theNearest3;
00108 int theHardest;
00109
00110 Handle <GenParticleRefVector> particles;
00111
00112 edm::InputTag m_jetsSrc, m_ParticleSrc;
00113 double coneSizeToAssociate;
00114 bool physDefinition;
00115 bool doPriority;
00116 vector<int> priorityList;
00117
00118 };
00119
00120
00121
00122 JetPartonMatcher::JetPartonMatcher( const edm::ParameterSet& iConfig ) :
00123 theHeaviest(0),
00124 theNearest2(0),
00125 theNearest3(0),
00126 theHardest(0)
00127 {
00128 produces<JetMatchedPartonsCollection>();
00129 m_jetsSrc = iConfig.getParameter<edm::InputTag>("jets");
00130 m_ParticleSrc = iConfig.getParameter<edm::InputTag>("partons");
00131 coneSizeToAssociate = iConfig.getParameter<double>("coneSizeToAssociate");
00132 if ( iConfig.exists("doPriority") ) {
00133 doPriority = iConfig.getParameter<bool>("doPriority");
00134 priorityList = iConfig.getParameter<vector<int> >("priorityList");
00135 } else {
00136 doPriority = false;
00137 priorityList.clear();
00138 }
00139 }
00140
00141
00142
00143 JetPartonMatcher::~JetPartonMatcher()
00144 {
00145 }
00146
00147
00148
00149 void JetPartonMatcher::produce( Event& iEvent, const EventSetup& iEs )
00150 {
00151 edm::Handle <edm::View <reco::Jet> > jets_h;
00152 iEvent.getByLabel(m_jetsSrc, jets_h );
00153 iEvent.getByLabel(m_ParticleSrc, particles );
00154
00155 edm::LogVerbatim("JetPartonMatcher") << "=== Partons size:" << particles->size();
00156
00157 for( size_t m = 0; m != particles->size(); ++ m ) {
00158 const GenParticle & aParton = *(particles->at(m).get());
00159 edm::LogVerbatim("JetPartonMatcher") << aParton.status() << " " <<
00160 aParton.pdgId() << " " <<
00161 aParton.pt() << " " <<
00162 aParton.eta() << " " <<
00163 aParton.phi() << endl;
00164 }
00165
00166 auto_ptr<reco::JetMatchedPartonsCollection> jetMatchedPartons( new JetMatchedPartonsCollection(reco::JetRefBaseProd(jets_h)));
00167
00168 for (size_t j = 0; j < jets_h->size(); j++) {
00169
00170 const int theMappedPartonAlg = fillAlgoritDefinition( (*jets_h)[j] );
00171 const int theMappedPartonPhy = fillPhysicsDefinition( (*jets_h)[j] );
00172
00173 GenParticleRef pHV;
00174 GenParticleRef pN2;
00175 GenParticleRef pN3;
00176 GenParticleRef pPH;
00177 GenParticleRef pAL;
00178
00179 if(theHeaviest>=0) pHV = particles->at( theHeaviest );
00180 if(theNearest2>=0) pN2 = particles->at( theNearest2 );
00181 if(theNearest3>=0) pN3 = particles->at( theNearest3 );
00182 if(theMappedPartonPhy>=0) pPH = particles->at( theMappedPartonPhy );
00183 if(theMappedPartonAlg>=0) pAL = particles->at( theMappedPartonAlg );
00184
00185 (*jetMatchedPartons)[jets_h->refAt(j)]=MatchedPartons(pHV,pN2,pN3,pPH,pAL);
00186 }
00187
00188 iEvent.put( jetMatchedPartons );
00189
00190 }
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215 int JetPartonMatcher::fillAlgoritDefinition( const Jet& theJet ) {
00216
00217 int tempParticle = -1;
00218 int tempPartonHighestPt = -1;
00219 int tempNearest = -1;
00220 float maxPt = 0;
00221 float minDr = 1000;
00222 bool foundPriority = false;
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243 for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
00244 const Candidate & aParton = *(particles->at(m).get());
00245
00246
00247
00248
00249 if ( doPriority ) {
00250 int ipdgId = aParton.pdgId();
00251 vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
00252
00253 if ( ipriority != priorityList.end() ) {
00254
00255
00256 double dist = DeltaR( theJet.p4(), aParton.p4() );
00257 if( dist <= coneSizeToAssociate ) {
00258 tempParticle = m;
00259 foundPriority = true;
00260 }
00261 }
00262 }
00263
00264
00265 else {
00266 foundPriority = false;
00267 }
00268
00269
00270
00271
00272
00273
00274 if( !foundPriority && aParton.status() != 3 ) {
00275 double dist = DeltaR( theJet.p4(), aParton.p4() );
00276 if( dist <= coneSizeToAssociate ) {
00277 if( dist < minDr ) {
00278 minDr = dist;
00279 tempNearest = m;
00280 }
00281 if( tempParticle == -1 && ( abs( aParton.pdgId() ) == 4 ) ) tempParticle = m;
00282 if( abs( aParton.pdgId() ) == 5 ) tempParticle = m;
00283 if( aParton.pt() > maxPt ) {
00284 maxPt = aParton.pt();
00285 tempPartonHighestPt = m;
00286 }
00287 }
00288 }
00289 }
00290
00291 if ( foundPriority ) {
00292 theHeaviest = tempParticle;
00293 theHardest = -1;
00294 theNearest2 = -1;
00295 } else {
00296 theHeaviest = tempParticle;
00297 theHardest = tempPartonHighestPt;
00298 theNearest2 = tempNearest;
00299 if ( tempParticle == -1 ) tempParticle = tempPartonHighestPt;
00300 }
00301 return tempParticle;
00302 }
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325
00326 int JetPartonMatcher::fillPhysicsDefinition( const Jet& theJet ) {
00327
00328 float TheBiggerConeSize = 0.7;
00329 int tempParticle = -1;
00330 int nInTheCone = 0;
00331 int tempNearest = -1;
00332 float minDr = 1000;
00333 bool foundPriority = false;
00334
00335 vector<const reco::Candidate *> theContaminations;
00336 theContaminations.clear();
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353
00354
00355
00356
00357 for( size_t m = 0; m != particles->size() && !foundPriority; ++ m ) {
00358
00359 const Candidate & aParticle = *(particles->at(m).get());
00360
00361
00362
00363
00364 if ( doPriority ) {
00365 int ipdgId = aParticle.pdgId();
00366 vector<int>::const_iterator ipriority = find( priorityList.begin(), priorityList.end(), abs(ipdgId) );
00367
00368 if ( ipriority != priorityList.end() ) {
00369
00370
00371 double dist = DeltaR( theJet.p4(), aParticle.p4() );
00372 if( dist <= coneSizeToAssociate ) {
00373 tempParticle = m;
00374 foundPriority = true;
00375 }
00376 }
00377 }
00378
00379
00380 else {
00381 foundPriority = false;
00382 }
00383
00384
00385 if ( !foundPriority ) {
00386
00387
00388 bool isAParton = false;
00389 int flavour = abs(aParticle.pdgId());
00390 if(flavour == 1 ||
00391 flavour == 2 ||
00392 flavour == 3 ||
00393 flavour == 4 ||
00394 flavour == 5 ||
00395 flavour == 21 ) isAParton = true;
00396 if(!isAParton) continue;
00397 double dist = DeltaR( theJet.p4(), aParticle.p4() );
00398 if( aParticle.status() == 3 && dist < minDr ) {
00399 minDr = dist;
00400 tempNearest = m;
00401 }
00402 if( aParticle.status() == 3 && dist <= coneSizeToAssociate ) {
00403
00404 tempParticle = m;
00405 nInTheCone++;
00406 }
00407
00408 if( aParticle.numberOfDaughters() > 0 && aParticle.status() != 3 ) {
00409 if( flavour == 1 ||
00410 flavour == 2 ||
00411 flavour == 3 ||
00412 flavour == 21 ) continue;
00413 if( dist < TheBiggerConeSize ) theContaminations.push_back( &aParticle );
00414 }
00415 }
00416 }
00417
00418
00419
00420 if ( !foundPriority ) {
00421 theNearest3 = tempNearest;
00422
00423 if(nInTheCone != 1) return -1;
00424 if(theContaminations.size() == 0 ) return tempParticle;
00425 int initialPartonFlavour = abs( (particles->at(tempParticle).get()) ->pdgId() );
00426
00427 vector<const Candidate *>::const_iterator itCont = theContaminations.begin();
00428 for( ; itCont != theContaminations.end(); itCont++ ) {
00429 int contaminatingFlavour = abs( (*itCont)->pdgId() );
00430 if( (*itCont)->numberOfMothers()>0 && (*itCont)->mother(0) == particles->at(tempParticle).get() ) continue;
00431 if( initialPartonFlavour == 4 ) {
00432 if( contaminatingFlavour == 4 ) continue;
00433 tempParticle = -1;
00434 return tempParticle;
00435 }
00436 }
00437 }
00438
00439 else {
00440 theHeaviest = tempParticle;
00441 theNearest2 = -1;
00442 theNearest3 = -1;
00443 theHardest = -1;
00444 }
00445
00446 return tempParticle;
00447 }
00448
00449
00450 DEFINE_FWK_MODULE(JetPartonMatcher);
00451