CMS 3D CMS Logo

MCMultiParticleFilter.cc
Go to the documentation of this file.
3 
5  src_(iConfig.getUntrackedParameter<edm::InputTag>("src",edm::InputTag(std::string("generator"),"unsmeared"))),
6  token_(consumes<edm::HepMCProduct>(src_)),
7  numRequired_(iConfig.getParameter<int>("NumRequired")),
8  acceptMore_(iConfig.getParameter<bool>("AcceptMore")),
9  particleID_(iConfig.getParameter< std::vector<int> >("ParticleID")),
10  ptMin_(iConfig.getParameter< std::vector<double> >("PtMin")),
11  etaMax_(iConfig.getParameter< std::vector<double> >("EtaMax")),
12  status_(iConfig.getParameter< std::vector<int> >("Status")),
13  totalEvents_(0), passedEvents_(0)
14 {
15  //here do whatever other initialization is needed
16 
17  // default pt, eta, status cuts to "don't care"
18  std::vector<double> defptmin(1, 0);
19  std::vector<double> defetamax(1, 999.0);
20  std::vector<int> defstat(1, 0);
21  std::vector<int> defmother;
22  defmother.push_back(0);
23  motherID_ = iConfig.getUntrackedParameter< std::vector<int> >("MotherID", defmother);
24 
25  std::vector<double> defDecayRadiusmin;
26  defDecayRadiusmin.push_back(-1.);
27  decayRadiusMin = iConfig.getUntrackedParameter<std::vector<double> >("MinDecayRadius", defDecayRadiusmin);
28 
29  std::vector<double> defDecayRadiusmax;
30  defDecayRadiusmax.push_back(1.e5);
31  decayRadiusMax = iConfig.getUntrackedParameter<std::vector<double> >("MaxDecayRadius", defDecayRadiusmax);
32 
33  std::vector<double> defDecayZmin;
34  defDecayZmin.push_back(-1.e5);
35  decayZMin = iConfig.getUntrackedParameter<std::vector<double> >("MinDecayZ", defDecayZmin);
36 
37  std::vector<double> defDecayZmax;
38  defDecayZmax.push_back(1.e5);
39  decayZMax = iConfig.getUntrackedParameter<std::vector<double> >("MaxDecayZ", defDecayZmax);
40 
41  // check for same size
42  if ( (ptMin_.size() > 1 && particleID_.size() != ptMin_.size())
43  || (etaMax_.size() > 1 && particleID_.size() != etaMax_.size())
44  || (status_.size() > 1 && particleID_.size() != status_.size())
45  || (motherID_.size() > 1 && particleID_.size() != motherID_.size())
46  || (decayRadiusMin.size() > 1 && particleID_.size() != decayRadiusMin.size())
47  || (decayRadiusMax.size() > 1 && particleID_.size() != decayRadiusMax.size())
48  || (decayZMin.size() > 1 && particleID_.size() != decayZMin.size())
49  || (decayZMax.size() > 1 && particleID_.size() != decayZMax.size())
50  ) {
51  edm::LogWarning("MCMultiParticleFilter") << "WARNING: MCMultiParticleFilter: size of PtMin, EtaMax, motherID, decayRadiusMin, decayRadiusMax, decayZMin, decayZMax and/or Status does not match ParticleID size!" << std::endl;
52  }
53 
54  // Fill arrays with defaults if necessary
55  while (ptMin_.size() < particleID_.size())
56  ptMin_.push_back(defptmin[0]);
57  while (etaMax_.size() < particleID_.size())
58  etaMax_.push_back(defetamax[0]);
59  while (status_.size() < particleID_.size())
60  status_.push_back(defstat[0]);
61  while (motherID_.size() < particleID_.size())
62  motherID_.push_back(defmother[0]);
63 
64  // if decayRadiusMin size smaller than particleID , fill up further with defaults
65  if (particleID_.size() > decayRadiusMin.size()) {
66  for (unsigned int i = decayRadiusMin.size(); i < particleID_.size(); i++) {
67  decayRadiusMin.push_back(-10.);
68  }
69  }
70  // if decayRadiusMax size smaller than particleID , fill up further with defaults
71  if (particleID_.size() > decayRadiusMax.size()) {
72  for (unsigned int i = decayRadiusMax.size(); i < particleID_.size(); i++) {
73  decayRadiusMax.push_back(1.e5);
74  }
75  }
76 
77  // if decayZMin size smaller than particleID , fill up further with defaults
78  if (particleID_.size() > decayZMin.size()) {
79  for (unsigned int i = decayZMin.size(); i < particleID_.size(); i++) {
80  decayZMin.push_back(-1.e5);
81  }
82  }
83  // if decayZMax size smaller than particleID , fill up further with defaults
84  if (particleID_.size() > decayZMax.size()) {
85  for (unsigned int i = decayZMax.size(); i < particleID_.size(); i++) {
86  decayZMax.push_back(1.e5);
87  }
88  }
89 }
90 
92 {
93 
94  // do anything here that needs to be done at destruction time
95  // (e.g. close files, deallocate resources etc.)
96 
97 }
98 
99 
100 // ------------ method called to skim the data ------------
102 {
104  iEvent.getByToken(token_, evt);
105 
106  totalEvents_++;
107  int nFound = 0;
108 
109  const HepMC::GenEvent * myGenEvent = evt->GetEvent();
110 
111  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
112  p != myGenEvent->particles_end(); ++p ) {
113 
114  for (unsigned int i = 0; i < particleID_.size(); ++i) {
115  if ((particleID_[i] == 0 || abs(particleID_[i]) == abs((*p)->pdg_id())) &&
116  (*p)->momentum().perp() > ptMin_[i] &&
117  fabs((*p)->momentum().eta()) < etaMax_[i] &&
118  (status_[i] == 0 || (*p)->status() == status_[i])) {
119  if (!((*p)->production_vertex()))
120  continue;
121 
122  double decx = (*p)->production_vertex()->position().x();
123  double decy = (*p)->production_vertex()->position().y();
124  double decrad = sqrt(decx * decx + decy * decy);
125  if (decrad < decayRadiusMin[i])
126  continue;
127  if (decrad > decayRadiusMax[i])
128  continue;
129 
130  double decz = (*p)->production_vertex()->position().z();
131  if (decz < decayZMin[i])
132  continue;
133  if (decz > decayZMax[i])
134  continue;
135 
136  if(motherID_[i] == 0 ){ // do not check for mother ID if not sepcified
137  nFound++;
138  break; // only match a given particle once!
139  }
140  else{
141  bool hascorrectmother=false;
142  for ( HepMC::GenVertex::particles_in_const_iterator mo = (*p)->production_vertex()->particles_in_const_begin(); mo != (*p)->production_vertex()->particles_in_const_end(); ++mo){
143  if( (*mo)->pdg_id() == motherID_[i]){
144  hascorrectmother = true;
145  break;
146  }
147  }
148  if(hascorrectmother){
149  nFound++;
150  break; // only match a given particle once!
151  }
152  }
153  }
154  } // loop over targets
155 
156  if (acceptMore_ && nFound == numRequired_) break; // stop looking if we don't mind having more
157  } // loop over particles
158 
159  if (nFound == numRequired_) {
160  passedEvents_++;
161  return true;
162  } else {
163  return false;
164  }
165 
166 }
167 
168 // ------------ method called once each job just after ending the event loop ------------
170  edm::LogInfo("MCMultiParticleFilter") << "=== Results of MCMultiParticleFilter: passed "
171  << passedEvents_ << "/" << totalEvents_ << " events" << std::endl;
172 }
173 
174 //define this as a plug-in
176 
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< edm::HepMCProduct > token_
bool filter(edm::Event &, const edm::EventSetup &) override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< double > decayRadiusMin
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< int > status_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > decayZMax
std::vector< int > particleID_
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
std::vector< double > decayRadiusMax
HLT enums.
std::vector< int > motherID_
MCMultiParticleFilter(const edm::ParameterSet &)
std::vector< double > decayZMin
std::vector< double > ptMin_
std::vector< double > etaMax_