CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LHEGenericFilter.cc
Go to the documentation of this file.
2 
3 using namespace edm;
4 using namespace std;
5 
7 numRequired_(iConfig.getParameter<int>("NumRequired")),
8 acceptLogic_(iConfig.getParameter<std::string>("AcceptLogic")),
9 particleID_(iConfig.getParameter< std::vector<int> >("ParticleID")),
10 totalEvents_(0), passedEvents_(0)
11 {
12  //here do whatever other initialization is needed
13  src_ = consumes<LHEEventProduct>(iConfig.getParameter<edm::InputTag>("src"));
14 
15  if(acceptLogic_.compare("LT")==0) whichlogic = LT;
16  else if(acceptLogic_.compare("GT")==0) whichlogic = GT;
17  else if(acceptLogic_.compare("EQ")==0) whichlogic = EQ;
18  else if(acceptLogic_.compare("NE")==0) whichlogic = NE;
19  else edm::LogError ("cat_A") << "wrong input for AcceptLogic string";
20 
21 
22 }
23 
25 {
26 
27  // do anything here that needs to be done at destruction time
28  // (e.g. close files, deallocate resources etc.)
29 
30 }
31 
32 
33 // ------------ method called to skim the data ------------
35 {
37  iEvent.getByToken( src_ , EvtHandle ) ;
38 
39  totalEvents_++;
40  int nFound = 0;
41 
42  for (int i = 0; i < EvtHandle->hepeup().NUP; ++i) {
43  if (EvtHandle->hepeup().ISTUP[i] != 1) { // keep only outgoing particles
44  continue;
45  }
46  for (unsigned int j = 0; j < particleID_.size(); ++j) {
47  if (particleID_[j] == 0 || abs(particleID_[j]) == abs(EvtHandle->hepeup().IDUP[i]) ) {
48  nFound++;
49  break; // only match a given particle once!
50  }
51  } // loop over targets
52 
53  } // loop over particles
54 
55  // event accept/reject logic
56  if (
57  (whichlogic==LT && nFound < numRequired_)
58  || (whichlogic==GT && nFound > numRequired_)
59  || (whichlogic==EQ && nFound == numRequired_)
60  || (whichlogic==NE && nFound != numRequired_)
61  ) {
62  passedEvents_++;
63  return true;
64  } else {
65  return false;
66  }
67 
68 }
69 
70 // ------------ method called once each job just after ending the event loop ------------
72  edm::LogInfo("LHEGenericFilter") << "=== Results of LHEGenericFilter: passed "
73  << passedEvents_ << "/" << totalEvents_ << " events" << std::endl;
74 }
75 
76 //define this as a plug-in
78 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int iEvent
Definition: GenABIO.cc:230
virtual bool filter(edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< LHEEventProduct > src_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
virtual void endJob()
std::string acceptLogic_
LHEGenericFilter(const edm::ParameterSet &)
std::vector< int > particleID_