Skip to content

Commit

Permalink
Merge pull request #11149 from perrozzi/from-CMSSW_7_1_19
Browse files Browse the repository at this point in the history
Changes to PartonShowerBsHepMCFilter and LHEGenericFilter
  • Loading branch information
davidlange6 committed Sep 22, 2015
2 parents 3b94c4e + ddea33b commit bdc86d3
Show file tree
Hide file tree
Showing 7 changed files with 128 additions and 17 deletions.
1 change: 1 addition & 0 deletions GeneratorInterface/Core/BuildFile.xml
Original file line number Diff line number Diff line change
Expand Up @@ -3,6 +3,7 @@
<use name="FWCore/Utilities"/>
<use name="SimDataFormats/GeneratorProducts"/>
<use name="GeneratorInterface/LHEInterface"/>
<use name="heppdt"/>
<use name="boost"/>
<use name="clhep"/>
<use name="lhapdf"/>
Expand Down
30 changes: 30 additions & 0 deletions GeneratorInterface/Core/interface/PartonShowerBsHepMCFilter.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,30 @@
/**
** Description: Filter gen particles based on pdg_id and status code
**
** @author bortigno
** @version 1.0 02.04.2015
*/


#ifndef __PARTONSHOWERBSHEPMCFILTER__
#define __PARTONSHOWERBSHEPMCFILTER__


#include "GeneratorInterface/Core/interface/BaseHepMCFilter.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"

class PartonShowerBsHepMCFilter : public BaseHepMCFilter{

public:

PartonShowerBsHepMCFilter( const edm::ParameterSet & );
~PartonShowerBsHepMCFilter();

virtual bool filter(const HepMC::GenEvent* evt);

private:

};


#endif
5 changes: 5 additions & 0 deletions GeneratorInterface/Core/src/HepMCFilterDriver.cc
Original file line number Diff line number Diff line change
@@ -1,5 +1,7 @@
#include "GeneratorInterface/Core/interface/HepMCFilterDriver.h"
#include "GeneratorInterface/Core/interface/GenericDauHepMCFilter.h"
#include "GeneratorInterface/Core/interface/PartonShowerBsHepMCFilter.h"


#include "FWCore/MessageLogger/interface/MessageLogger.h"

Expand All @@ -21,6 +23,9 @@ HepMCFilterDriver::HepMCFilterDriver(const edm::ParameterSet& pset) :
if (filterName=="GenericDauHepMCFilter") {
filter_ = new GenericDauHepMCFilter(filterParameters);
}
else if (filterName=="PartonShowerBsHepMCFilter") {
filter_ = new PartonShowerBsHepMCFilter(filterParameters);
}
else {
edm::LogError("HepMCFilterDriver")<< "Invalid HepMCFilter name:" << filterName;
}
Expand Down
48 changes: 48 additions & 0 deletions GeneratorInterface/Core/src/PartonShowerBsHepMCFilter.cc
Original file line number Diff line number Diff line change
@@ -0,0 +1,48 @@
#include "GeneratorInterface/Core/interface/PartonShowerBsHepMCFilter.h"
#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
#include <iostream>
#include "HepPDT/ParticleID.hh"


using namespace edm;
using namespace std;


//constructor
PartonShowerBsHepMCFilter::PartonShowerBsHepMCFilter(const edm::ParameterSet& iConfig)
{

}


//destructor
PartonShowerBsHepMCFilter::~PartonShowerBsHepMCFilter()
{

}

//
// member functions
//

// ------------ method called to produce the data ------------
bool PartonShowerBsHepMCFilter::filter(const HepMC::GenEvent* evt)
{

// loop over gen particles
for ( HepMC::GenEvent::particle_const_iterator p = evt->particles_begin(); p != evt->particles_end(); ++p ){

// check only status 2 particles
if( (*p)->status()==2 ){
// if one of the status 2 particles is a B-hadron, accept the event
HepPDT::ParticleID pid((*p)->pdg_id());
if( pid.hasBottom() ){
return true; // accept event
}
}

}

return false; // skip event

}
8 changes: 6 additions & 2 deletions GeneratorInterface/GenFilters/interface/LHEGenericFilter.h
Original file line number Diff line number Diff line change
Expand Up @@ -49,10 +49,14 @@ class LHEGenericFilter : public edm::EDFilter {

edm::EDGetTokenT<LHEEventProduct> src_;
int numRequired_; // number of particles required to pass filter
bool acceptMore_; // if true (default), accept numRequired or more.
// if false, accept events with exactly equal to numRequired.
std::string acceptLogic_; // LT meaning <
// GT >
// EQ =
// NE !=
std::vector<int> particleID_; // vector of particle IDs to look for
int totalEvents_; // counters
int passedEvents_;
enum logic_ { LT, GT, EQ, NE};
logic_ whichlogic;
};
#endif
8 changes: 8 additions & 0 deletions GeneratorInterface/GenFilters/python/LHEGenericFilter_cfi.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,8 @@
import FWCore.ParameterSet.Config as cms

lheGenericFilter = cms.EDFilter("LHEGenericFilter",
src = cms.InputTag("source"),
NumRequired = cms.int32(2),
ParticleID = cms.vint32(5),
AcceptLogic = cms.string("LT") # LT meaning < NumRequired, GT >, EQ =, NE !=
)
45 changes: 30 additions & 15 deletions GeneratorInterface/GenFilters/src/LHEGenericFilter.cc
Original file line number Diff line number Diff line change
@@ -1,20 +1,30 @@
#include "GeneratorInterface/GenFilters/interface/LHEGenericFilter.h"

using namespace edm;
using namespace std;

LHEGenericFilter::LHEGenericFilter(const edm::ParameterSet& iConfig) :
numRequired_(iConfig.getParameter<int>("NumRequired")),
acceptMore_(iConfig.getParameter<bool>("AcceptMore")),
particleID_(iConfig.getParameter< std::vector<int> >("ParticleID")),
totalEvents_(0), passedEvents_(0)
numRequired_(iConfig.getParameter<int>("NumRequired")),
acceptLogic_(iConfig.getParameter<std::string>("AcceptLogic")),
particleID_(iConfig.getParameter< std::vector<int> >("ParticleID")),
totalEvents_(0), passedEvents_(0)
{
//here do whatever other initialization is needed
src_ = consumes<LHEEventProduct>(iConfig.getParameter<edm::InputTag>("src"));

if(acceptLogic_.compare("LT")==0) whichlogic = LT;
else if(acceptLogic_.compare("GT")==0) whichlogic = GT;
else if(acceptLogic_.compare("EQ")==0) whichlogic = EQ;
else if(acceptLogic_.compare("NE")==0) whichlogic = NE;
else edm::LogError ("cat_A") << "wrong input for AcceptLogic string";

}

LHEGenericFilter::~LHEGenericFilter()
{
// do anything here that needs to be done at destruction time
// (e.g. close files, deallocate resources etc.)

// do anything here that needs to be done at destruction time
// (e.g. close files, deallocate resources etc.)

}

Expand All @@ -29,20 +39,25 @@ bool LHEGenericFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
int nFound = 0;

for (int i = 0; i < EvtHandle->hepeup().NUP; ++i) {
if (EvtHandle->hepeup().ISTUP[i] != 1) {
continue;
}
if (EvtHandle->hepeup().ISTUP[i] != 1) { // keep only outgoing particles
continue;
}
for (unsigned int j = 0; j < particleID_.size(); ++j) {
if (particleID_[j] == 0 || abs(particleID_[j]) == abs(EvtHandle->hepeup().IDUP[i]) ) {
nFound++;
break; // only match a given particle once!
nFound++;
break; // only match a given particle once!
}
} // loop over targets

if (acceptMore_ && nFound == numRequired_) break; // stop looking if we don't mind having more
} // loop over particles

if (nFound == numRequired_) {
// event accept/reject logic
if (
(whichlogic==LT && nFound < numRequired_)
|| (whichlogic==GT && nFound > numRequired_)
|| (whichlogic==EQ && nFound == numRequired_)
|| (whichlogic==NE && nFound != numRequired_)
) {
passedEvents_++;
return true;
} else {
Expand All @@ -54,7 +69,7 @@ bool LHEGenericFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
// ------------ method called once each job just after ending the event loop ------------
void LHEGenericFilter::endJob() {
edm::LogInfo("LHEGenericFilter") << "=== Results of LHEGenericFilter: passed "
<< passedEvents_ << "/" << totalEvents_ << " events" << std::endl;
<< passedEvents_ << "/" << totalEvents_ << " events" << std::endl;
}

//define this as a plug-in
Expand Down

0 comments on commit bdc86d3

Please sign in to comment.