86 #include <Math/VectorUtil.h> 127 produces<JetMatchedPartonsCollection>();
131 if ( iConfig.
exists(
"doPriority") ) {
133 priorityList = iConfig.
getParameter<vector<int> >(
"priorityList");
136 priorityList.clear();
160 aParton.
pdgId() <<
" " <<
161 aParton.
pt() <<
" " <<
162 aParton.
eta() <<
" " <<
163 aParton.
phi() << endl;
166 auto jetMatchedPartons = std::make_unique<JetMatchedPartonsCollection>(
reco::JetRefBaseProd(jets_h));
168 for (
size_t j = 0; j < jets_h->size(); j++) {
170 const int theMappedPartonAlg = fillAlgoritDefinition( (*jets_h)[j], wv );
171 const int theMappedPartonPhy = fillPhysicsDefinition( (*jets_h)[j], wv );
182 if(theMappedPartonPhy>=0) pPH = wv.
particles->at( theMappedPartonPhy );
183 if(theMappedPartonAlg>=0) pAL = wv.
particles->at( theMappedPartonAlg );
185 (*jetMatchedPartons)[jets_h->refAt(j)]=
MatchedPartons(pHV,pN2,pN3,pPH,pAL);
217 int tempParticle = -1;
218 int tempPartonHighestPt = -1;
219 int tempNearest = -1;
222 bool foundPriority =
false;
243 for(
size_t m = 0;
m != wv.
particles->size() && !foundPriority; ++
m ) {
250 int ipdgId = aParton.
pdgId();
251 vector<int>::const_iterator ipriority =
find( priorityList.begin(), priorityList.end(),
abs(ipdgId) );
253 if ( ipriority != priorityList.end() ) {
256 double dist =
DeltaR( theJet.
p4(), aParton.
p4() );
259 foundPriority =
true;
266 foundPriority =
false;
274 if( !foundPriority && aParton.
status() != 3 ) {
275 double dist =
DeltaR( theJet.
p4(), aParton.
p4() );
281 if( tempParticle == -1 && (
abs( aParton.
pdgId() ) == 4 ) ) tempParticle =
m;
282 if(
abs( aParton.
pdgId() ) == 5 ) tempParticle =
m;
284 maxPt = aParton.
pt();
285 tempPartonHighestPt =
m;
291 if ( foundPriority ) {
299 if ( tempParticle == -1 ) tempParticle = tempPartonHighestPt;
328 float TheBiggerConeSize = 0.7;
329 int tempParticle = -1;
331 int tempNearest = -1;
333 bool foundPriority =
false;
335 vector<const reco::Candidate *> theContaminations;
336 theContaminations.clear();
357 for(
size_t m = 0;
m != wv.
particles->size() && !foundPriority; ++
m ) {
365 int ipdgId = aParticle.
pdgId();
366 vector<int>::const_iterator ipriority =
find( priorityList.begin(), priorityList.end(),
abs(ipdgId) );
368 if ( ipriority != priorityList.end() ) {
371 double dist =
DeltaR( theJet.
p4(), aParticle.
p4() );
374 foundPriority =
true;
381 foundPriority =
false;
385 if ( !foundPriority ) {
388 bool isAParton =
false;
395 flavour == 21 ) isAParton =
true;
396 if(!isAParton)
continue;
397 double dist =
DeltaR( theJet.
p4(), aParticle.
p4() );
398 if( aParticle.
status() == 3 && dist < minDr ) {
412 flavour == 21 )
continue;
413 if( dist < TheBiggerConeSize ) theContaminations.push_back( &aParticle );
420 if ( !foundPriority ) {
423 if(nInTheCone != 1)
return -1;
424 if(theContaminations.size() == 0 )
return tempParticle;
425 int initialPartonFlavour =
abs( (wv.
particles->at(tempParticle).get()) ->
pdgId() );
427 vector<const Candidate *>::const_iterator itCont = theContaminations.begin();
428 for( ; itCont != theContaminations.end(); itCont++ ) {
429 int contaminatingFlavour =
abs( (*itCont)->pdgId() );
430 if( (*itCont)->numberOfMothers()>0 && (*itCont)->mother(0) == wv.
particles->at(tempParticle).get() )
continue;
431 if( initialPartonFlavour == 4 ) {
432 if( contaminatingFlavour == 4 )
continue;
T getParameter(std::string const &) const
virtual double pt() const final
transverse momentum
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
virtual double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
#define DEFINE_FWK_MODULE(type)
int fillAlgoritDefinition(const Jet &, WorkingVariables &) const
bool exists(std::string const ¶meterName) const
checks if a parameter exists
virtual int status() const final
status word
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
virtual double phi() const final
momentum azimuthal angle
edm::RefToBaseProd< reco::Jet > JetRefBaseProd
virtual int status() const =0
status word
edm::EDGetTokenT< GenParticleRefVector > m_ParticleSrcToken
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
virtual int pdgId() const =0
PDG identifier.
virtual int pdgId() const final
PDG identifier.
Abs< T >::type abs(const T &t)
JetPartonMatcher(const edm::ParameterSet &)
vector< int > priorityList
virtual double pt() const =0
transverse momentum
double coneSizeToAssociate
int fillPhysicsDefinition(const Jet &, WorkingVariables &) const
edm::EDGetTokenT< edm::View< reco::Jet > > m_jetsSrcToken
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
virtual size_type numberOfDaughters() const =0
number of daughters
Handle< GenParticleRefVector > particles