CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MCMultiParticleFilter.cc
Go to the documentation of this file.
3 
5  src_(iConfig.getParameter<edm::InputTag>("src")),
6  numRequired_(iConfig.getParameter<int>("NumRequired")),
7  acceptMore_(iConfig.getParameter<bool>("AcceptMore")),
8  particleID_(iConfig.getParameter< std::vector<int> >("ParticleID")),
9  ptMin_(iConfig.getParameter< std::vector<double> >("PtMin")),
10  etaMax_(iConfig.getParameter< std::vector<double> >("EtaMax")),
11  status_(iConfig.getParameter< std::vector<int> >("Status")),
12  totalEvents_(0), passedEvents_(0)
13 {
14  //here do whatever other initialization is needed
15 
16  // default pt, eta, status cuts to "don't care"
17  std::vector<double> defptmin(1, 0);
18  std::vector<double> defetamax(1, 999.0);
19  std::vector<int> defstat(1, 0);
20 
21  // check for same size
22  if ( (ptMin_.size() > 1 && particleID_.size() != ptMin_.size())
23  || (etaMax_.size() > 1 && particleID_.size() != etaMax_.size())
24  || (status_.size() > 1 && particleID_.size() != status_.size()) ) {
25  edm::LogWarning("MCMultiParticleFilter") << "WARNING: MCMultiParticleFilter: size of PtMin, EtaMax, and/or Status does not match ParticleID size!" << std::endl;
26  }
27 
28  // Fill arrays with defaults if necessary
29  while (ptMin_.size() < particleID_.size())
30  ptMin_.push_back(defptmin[0]);
31  while (etaMax_.size() < particleID_.size())
32  etaMax_.push_back(defetamax[0]);
33  while (status_.size() < particleID_.size())
34  status_.push_back(defstat[0]);
35 }
36 
38 {
39 
40  // do anything here that needs to be done at destruction time
41  // (e.g. close files, deallocate resources etc.)
42 
43 }
44 
45 
46 // ------------ method called to skim the data ------------
48 {
50  iEvent.getByLabel(src_, evt);
51 
52  totalEvents_++;
53  int nFound = 0;
54 
55  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
56 
57  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
58  p != myGenEvent->particles_end(); ++p ) {
59 
60  for (unsigned int i = 0; i < particleID_.size(); ++i) {
61  if ((particleID_[i] == 0 || abs(particleID_[i]) == abs((*p)->pdg_id())) &&
62  (*p)->momentum().perp() > ptMin_[i] &&
63  fabs((*p)->momentum().eta()) < etaMax_[i] &&
64  (status_[i] == 0 || (*p)->status() == status_[i])) {
65  nFound++;
66  break; // only match a given particle once!
67  }
68  } // loop over targets
69 
70  if (acceptMore_ && nFound == numRequired_) break; // stop looking if we don't mind having more
71  } // loop over particles
72 
73  if (nFound == numRequired_) {
74  passedEvents_++;
75  return true;
76  } else {
77  return false;
78  }
79 
80 }
81 
82 // ------------ method called once each job just after ending the event loop ------------
84  edm::LogInfo("MCMultiParticleFilter") << "=== Results of MCMultiParticleFilter: passed "
85  << passedEvents_ << "/" << totalEvents_ << " events" << std::endl;
86 }
87 
88 //define this as a plug-in
90 
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
int iEvent
Definition: GenABIO.cc:230
virtual bool filter(edm::Event &, const edm::EventSetup &)
std::vector< int > status_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:420
std::vector< int > particleID_
MCMultiParticleFilter(const edm::ParameterSet &)
std::vector< double > ptMin_
std::vector< double > etaMax_