CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GenParticlePruner.cc
Go to the documentation of this file.
10 
11 namespace helper {
12  struct SelectCode {
13  enum KeepOrDrop { kKeep, kDrop };
14  enum FlagDepth { kNone, kFirst, kAll };
17  bool all_;
18  };
19 }
20 
22 public:
24 private:
25  void produce(edm::Event&, const edm::EventSetup&) override;
29  std::vector<std::string> selection_;
30  std::vector<std::pair<StringCutObjectSelector<reco::GenParticle>, helper::SelectCode> > select_;
31  std::vector<int> flags_;
32  std::vector<size_t> indices_;
33  void parse(const std::string & selection, helper::SelectCode & code, std::string & cut) const;
34  void flagDaughters(const reco::GenParticle &, int);
35  void flagMothers(const reco::GenParticle &, int);
36  void recursiveFlagDaughters(size_t, const reco::GenParticleCollection &, int, std::vector<size_t> &);
37  void recursiveFlagMothers(size_t, const reco::GenParticleCollection &, int, std::vector<size_t> &);
38  void addDaughterRefs(std::vector<size_t> &, reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
39  void addMotherRefs(std::vector<size_t> &, reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
40 };
41 
42 using namespace edm;
43 using namespace std;
44 using namespace reco;
45 
46 const int keep = 1, drop = -1;
47 
49  using namespace ::helper;
50  size_t f = selection.find_first_not_of(' ');
51  size_t n = selection.size();
52  string command;
53  char c;
54  for(; (c = selection[f]) != ' ' && f < n; ++f) {
55  command.push_back(c);
56  }
57  if(command[0] == '+') {
58  command.erase(0, 1);
59  if(command[0] == '+') {
60  command.erase(0, 1);
61  code.mothersDepth_ = SelectCode::kAll;
62  } else {
63  code.mothersDepth_ = SelectCode::kFirst;
64  }
65  } else
67 
68  if(command[command.size() - 1] == '+') {
69  command.erase(command.size() - 1);
70  if(command[command.size()-1] == '+') {
71  command.erase(command.size() - 1);
72  code.daughtersDepth_ = SelectCode::kAll;
73  } else {
74  code.daughtersDepth_ = SelectCode::kFirst;
75  }
76  } else
78 
79  if(command == "keep") code.keepOrDrop_ = SelectCode::kKeep;
80  else if(command == "drop") code.keepOrDrop_ = SelectCode::kDrop;
81  else {
83  << "invalid selection command: " << command << "\n" << endl;
84  }
85  for(; f < n; ++f) {
86  if(selection[f] != ' ') break;
87  }
88  cut = string(selection, f);
89  if(cut[0] == '*')
90  cut = string(cut, 0, cut.find_first_of(' '));
91  code.all_ = cut == "*";
92 }
93 
95  firstEvent_(true),
96  srcToken_(consumes<GenParticleCollection>(cfg.getParameter<InputTag>("src"))), keepOrDropAll_(drop),
97  selection_(cfg.getParameter<vector<string> >("select")) {
98  using namespace ::helper;
99  produces<GenParticleCollection>();
100 }
101 
103  GenParticleRefVector daughters = gen.daughterRefVector();
104  for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i)
105  flags_[i->key()] = keepOrDrop;
106 }
107 
108 void GenParticlePruner::flagMothers(const reco::GenParticle & gen, int keepOrDrop) {
109  GenParticleRefVector mothers = gen.motherRefVector();
110  for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i)
111  flags_[i->key()] = keepOrDrop;
112 }
113 
115  std::vector<size_t> & allIndices ) {
116  GenParticleRefVector daughters = src[index].daughterRefVector();
117  // avoid infinite recursion if the daughters are set to "this" particle.
118  size_t cachedIndex = index;
119  for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i) {
120  index = i->key();
121  // To also avoid infinite recursion if a "loop" is found in the daughter list,
122  // check to make sure the index hasn't already been added.
123  if ( find( allIndices.begin(), allIndices.end(), index ) == allIndices.end() ) {
124  allIndices.push_back( index );
125  if ( cachedIndex != index ) {
126  flags_[index] = keepOrDrop;
127  recursiveFlagDaughters(index, src, keepOrDrop, allIndices);
128  }
129  }
130  }
131 }
132 
134  std::vector<size_t> & allIndices ) {
135  GenParticleRefVector mothers = src[index].motherRefVector();
136  // avoid infinite recursion if the mothers are set to "this" particle.
137  size_t cachedIndex = index;
138  for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i) {
139  index = i->key();
140  // To also avoid infinite recursion if a "loop" is found in the daughter list,
141  // check to make sure the index hasn't already been added.
142  if ( find( allIndices.begin(), allIndices.end(), index ) == allIndices.end() ) {
143  allIndices.push_back( index );
144  if ( cachedIndex != index ) {
145  flags_[index] = keepOrDrop;
146  recursiveFlagMothers(index, src, keepOrDrop, allIndices);
147  }
148  }
149  }
150 }
151 
153  if (firstEvent_) {
154  PdgEntryReplacer rep(es);
155  for(vector<string>::const_iterator i = selection_.begin(); i != selection_.end(); ++i) {
156  string cut;
158  parse(*i, code, cut);
159  if(code.all_) {
160  if(i != selection_.begin())
162  << "selections \"keep *\" and \"drop *\" can be used only as first options. Here used in position # "
163  << (i - selection_.begin()) + 1 << "\n" << endl;
164  switch(code.keepOrDrop_) {
165  case ::helper::SelectCode::kDrop :
166  keepOrDropAll_ = drop; break;
167  case ::helper::SelectCode::kKeep :
169  };
170  } else {
171  cut = rep.replace(cut);
172  select_.push_back(make_pair(StringCutObjectSelector<GenParticle>(cut), code));
173  }
174  }
175  firstEvent_ = false;
176  }
177 
178  using namespace ::helper;
180  evt.getByToken(srcToken_, src);
181  const size_t n = src->size();
182  flags_.clear();
183  flags_.resize(n, keepOrDropAll_);
184  for(size_t j = 0; j < select_.size(); ++j) {
185  const pair<StringCutObjectSelector<GenParticle>, SelectCode> & sel = select_[j];
186  SelectCode code = sel.second;
187  const StringCutObjectSelector<GenParticle> & cut = sel.first;
188  for(size_t i = 0; i < n; ++i) {
189  const GenParticle & p = (*src)[i];
190  if(cut(p)) {
191  int keepOrDrop = keep;
192  switch(code.keepOrDrop_) {
193  case SelectCode::kKeep:
194  keepOrDrop = keep; break;
195  case SelectCode::kDrop:
196  keepOrDrop = drop;
197  };
198  flags_[i] = keepOrDrop;
199  std::vector<size_t> allIndicesDa;
200  std::vector<size_t> allIndicesMo;
201  switch(code.daughtersDepth_) {
202  case SelectCode::kAll :
203  recursiveFlagDaughters(i, *src, keepOrDrop, allIndicesDa); break;
204  case SelectCode::kFirst :
205  flagDaughters(p, keepOrDrop); break;
206  case SelectCode::kNone:
207  ;
208  };
209  switch(code.mothersDepth_) {
210  case SelectCode::kAll :
211  recursiveFlagMothers(i, *src, keepOrDrop, allIndicesMo); break;
212  case SelectCode::kFirst :
213  flagMothers(p, keepOrDrop); break;
214  case SelectCode::kNone:
215  ;
216  };
217  }
218  }
219  }
220  indices_.clear();
221  int counter = 0;
222  for(size_t i = 0; i < n; ++i) {
223  if(flags_[i] == keep) {
224  indices_.push_back(i);
225  flags_[i] = counter++;
226  }
227  }
228 
229  auto_ptr<GenParticleCollection> out(new GenParticleCollection);
231  out->reserve(counter);
232  for(vector<size_t>::const_iterator i = indices_.begin(); i != indices_.end(); ++i) {
233  size_t index = *i;
234  const GenParticle & gen = (*src)[index];
235  const LeafCandidate & part = gen;
236  out->push_back(GenParticle(part));
237  GenParticle & newGen = out->back();
238  // The "daIndxs" and "moIndxs" keep a list of the keys for the mother/daughter
239  // parentage/descendency. In some cases, a circular referencing is encountered,
240  // which would result in an infinite loop. The list is checked to
241  // avoid this.
242  vector<size_t> daIndxs;
243  addDaughterRefs(daIndxs, newGen, outRef, gen.daughterRefVector());
244  vector<size_t> moIndxs;
245  addMotherRefs(moIndxs, newGen, outRef, gen.motherRefVector());
246  }
247 
248  evt.put(out);
249 }
250 
251 
252 void GenParticlePruner::addDaughterRefs(vector<size_t> & daIndxs,
253  GenParticle& newGen, GenParticleRefProd outRef,
254  const GenParticleRefVector& daughters) const {
255  for(GenParticleRefVector::const_iterator j = daughters.begin();
256  j != daughters.end(); ++j) {
257  GenParticleRef dau = *j;
258  if ( find(daIndxs.begin(), daIndxs.end(), dau.key()) == daIndxs.end() ) {
259  int idx = flags_[dau.key()];
260  daIndxs.push_back( dau.key() );
261  if(idx > 0) {
262  GenParticleRef newDau(outRef, static_cast<size_t>(idx));
263  newGen.addDaughter(newDau);
264  } else {
265  const GenParticleRefVector daus = dau->daughterRefVector();
266  if(daus.size()>0) {
267  addDaughterRefs(daIndxs, newGen, outRef, daus);
268  }
269  }
270  }
271  }
272 }
273 
274 
275 
276 void GenParticlePruner::addMotherRefs(vector<size_t> & moIndxs,
277  GenParticle& newGen, GenParticleRefProd outRef,
278  const GenParticleRefVector& mothers) const {
280  j != mothers.end(); ++j) {
281  GenParticleRef mom = *j;
282  if ( find(moIndxs.begin(), moIndxs.end(), mom.key()) == moIndxs.end() ) {
283  int idx = flags_[mom.key()];
284  moIndxs.push_back( mom.key() );
285  if(idx > 0) {
286  GenParticleRef newMom(outRef, static_cast<size_t>(idx));
287  newGen.addMother(newMom);
288  } else {
289  const GenParticleRefVector moms = mom->motherRefVector();
290  if(moms.size()>0)
291  addMotherRefs(moIndxs, newGen, outRef, moms);
292  }
293  }
294  }
295 }
296 
298 
std::vector< std::pair< StringCutObjectSelector< reco::GenParticle >, helper::SelectCode > > select_
void flagDaughters(const reco::GenParticle &, int)
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
int i
Definition: DBlmapReader.cc:9
string rep
Definition: cuy.py:1188
std::vector< int > flags_
void addMotherRefs(std::vector< size_t > &, reco::GenParticle &, reco::GenParticleRefProd, const reco::GenParticleRefVector &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void flagMothers(const reco::GenParticle &, int)
selection
main part
Definition: corrVsCorr.py:98
std::string replace(const std::string &) const
void recursiveFlagMothers(size_t, const reco::GenParticleCollection &, int, std::vector< size_t > &)
const daughters & daughterRefVector() const
references to daughtes
GenParticlePruner(const edm::ParameterSet &)
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void recursiveFlagDaughters(size_t, const reco::GenParticleCollection &, int, std::vector< size_t > &)
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
std::vector< size_t > indices_
const int keep
const mothers & motherRefVector() const
references to mothers
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
int j
Definition: DBlmapReader.cc:9
double f[11][100]
void parse(const std::string &selection, helper::SelectCode &code, std::string &cut) const
void addMother(const typename mothers::value_type &)
add a daughter via a reference
RefProd< PROD > getRefBeforePut()
Definition: Event.h:128
const int drop
tuple out
Definition: dbtoconf.py:99
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
part
Definition: HCALResponse.h:20
key_type key() const
Accessor for product key.
Definition: Ref.h:266
edm::EDGetTokenT< reco::GenParticleCollection > srcToken_
static std::atomic< unsigned int > counter
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< std::string > selection_
void addDaughterRefs(std::vector< size_t > &, reco::GenParticle &, reco::GenParticleRefProd, const reco::GenParticleRefVector &) const