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.
11 
12 namespace helper {
13  struct SelectCode {
14  enum KeepOrDrop { kKeep, kDrop };
15  enum FlagDepth { kNone, kFirst, kAll };
18  bool all_;
19  };
20 }
21 
23 public:
25 private:
26  void produce(edm::Event&, const edm::EventSetup&) override;
30  std::vector<std::string> selection_;
31  std::vector<std::pair<StringCutObjectSelector<reco::GenParticle>, helper::SelectCode> > select_;
32  std::vector<int> flags_;
33  std::vector<size_t> indices_;
34  void parse(const std::string & selection, helper::SelectCode & code, std::string & cut) const;
35  void flagDaughters(const reco::GenParticle &, int);
36  void flagMothers(const reco::GenParticle &, int);
37  void recursiveFlagDaughters(size_t, const reco::GenParticleCollection &, int, std::vector<size_t> &);
38  void recursiveFlagMothers(size_t, const reco::GenParticleCollection &, int, std::vector<size_t> &);
39  void addDaughterRefs(std::vector<size_t> &, reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
40  void addMotherRefs(std::vector<size_t> &, reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
41 };
42 
43 using namespace edm;
44 using namespace std;
45 using namespace reco;
46 
47 const int keep = 1, drop = -1;
48 
50  using namespace ::helper;
51  size_t f = selection.find_first_not_of(' ');
52  size_t n = selection.size();
53  string command;
54  char c;
55  for(; (c = selection[f]) != ' ' && f < n; ++f) {
56  command.push_back(c);
57  }
58  if(command[0] == '+') {
59  command.erase(0, 1);
60  if(command[0] == '+') {
61  command.erase(0, 1);
62  code.mothersDepth_ = SelectCode::kAll;
63  } else {
64  code.mothersDepth_ = SelectCode::kFirst;
65  }
66  } else
68 
69  if(command[command.size() - 1] == '+') {
70  command.erase(command.size() - 1);
71  if(command[command.size()-1] == '+') {
72  command.erase(command.size() - 1);
73  code.daughtersDepth_ = SelectCode::kAll;
74  } else {
75  code.daughtersDepth_ = SelectCode::kFirst;
76  }
77  } else
79 
80  if(command == "keep") code.keepOrDrop_ = SelectCode::kKeep;
81  else if(command == "drop") code.keepOrDrop_ = SelectCode::kDrop;
82  else {
84  << "invalid selection command: " << command << "\n" << endl;
85  }
86  for(; f < n; ++f) {
87  if(selection[f] != ' ') break;
88  }
89  cut = string(selection, f);
90  if(cut[0] == '*')
91  cut = string(cut, 0, cut.find_first_of(' '));
92  code.all_ = cut == "*";
93 }
94 
96  firstEvent_(true),
97  srcToken_(consumes<GenParticleCollection>(cfg.getParameter<InputTag>("src"))), keepOrDropAll_(drop),
98  selection_(cfg.getParameter<vector<string> >("select")) {
99  using namespace ::helper;
100  produces<GenParticleCollection>();
101  produces<edm::Association<reco::GenParticleCollection> >();
102 }
103 
105  GenParticleRefVector daughters = gen.daughterRefVector();
106  for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i)
107  flags_[i->key()] = keepOrDrop;
108 }
109 
110 void GenParticlePruner::flagMothers(const reco::GenParticle & gen, int keepOrDrop) {
111  GenParticleRefVector mothers = gen.motherRefVector();
112  for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i)
113  flags_[i->key()] = keepOrDrop;
114 }
115 
117  std::vector<size_t> & allIndices ) {
118  GenParticleRefVector daughters = src[index].daughterRefVector();
119  // avoid infinite recursion if the daughters are set to "this" particle.
120  size_t cachedIndex = index;
121  for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i) {
122  index = i->key();
123  // To also avoid infinite recursion if a "loop" is found in the daughter list,
124  // check to make sure the index hasn't already been added.
125  if ( find( allIndices.begin(), allIndices.end(), index ) == allIndices.end() ) {
126  allIndices.push_back( index );
127  if ( cachedIndex != index ) {
128  flags_[index] = keepOrDrop;
129  recursiveFlagDaughters(index, src, keepOrDrop, allIndices);
130  }
131  }
132  }
133 }
134 
136  std::vector<size_t> & allIndices ) {
137  GenParticleRefVector mothers = src[index].motherRefVector();
138  // avoid infinite recursion if the mothers are set to "this" particle.
139  size_t cachedIndex = index;
140  for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i) {
141  index = i->key();
142  // To also avoid infinite recursion if a "loop" is found in the daughter list,
143  // check to make sure the index hasn't already been added.
144  if ( find( allIndices.begin(), allIndices.end(), index ) == allIndices.end() ) {
145  allIndices.push_back( index );
146  if ( cachedIndex != index ) {
147  flags_[index] = keepOrDrop;
148  recursiveFlagMothers(index, src, keepOrDrop, allIndices);
149  }
150  }
151  }
152 }
153 
155  if (firstEvent_) {
156  PdgEntryReplacer rep(es);
157  for(vector<string>::const_iterator i = selection_.begin(); i != selection_.end(); ++i) {
158  string cut;
160  parse(*i, code, cut);
161  if(code.all_) {
162  if(i != selection_.begin())
164  << "selections \"keep *\" and \"drop *\" can be used only as first options. Here used in position # "
165  << (i - selection_.begin()) + 1 << "\n" << endl;
166  switch(code.keepOrDrop_) {
167  case ::helper::SelectCode::kDrop :
168  keepOrDropAll_ = drop; break;
169  case ::helper::SelectCode::kKeep :
171  };
172  } else {
173  cut = rep.replace(cut);
174  select_.push_back(make_pair(StringCutObjectSelector<GenParticle>(cut), code));
175  }
176  }
177  firstEvent_ = false;
178  }
179 
180  using namespace ::helper;
182  evt.getByToken(srcToken_, src);
183  const size_t n = src->size();
184  flags_.clear();
185  flags_.resize(n, keepOrDropAll_);
186  for(size_t j = 0; j < select_.size(); ++j) {
187  const pair<StringCutObjectSelector<GenParticle>, SelectCode> & sel = select_[j];
188  SelectCode code = sel.second;
189  const StringCutObjectSelector<GenParticle> & cut = sel.first;
190  for(size_t i = 0; i < n; ++i) {
191  const GenParticle & p = (*src)[i];
192  if(cut(p)) {
193  int keepOrDrop = keep;
194  switch(code.keepOrDrop_) {
195  case SelectCode::kKeep:
196  keepOrDrop = keep; break;
197  case SelectCode::kDrop:
198  keepOrDrop = drop;
199  };
200  flags_[i] = keepOrDrop;
201  std::vector<size_t> allIndicesDa;
202  std::vector<size_t> allIndicesMo;
203  switch(code.daughtersDepth_) {
204  case SelectCode::kAll :
205  recursiveFlagDaughters(i, *src, keepOrDrop, allIndicesDa); break;
206  case SelectCode::kFirst :
207  flagDaughters(p, keepOrDrop); break;
208  case SelectCode::kNone:
209  ;
210  };
211  switch(code.mothersDepth_) {
212  case SelectCode::kAll :
213  recursiveFlagMothers(i, *src, keepOrDrop, allIndicesMo); break;
214  case SelectCode::kFirst :
215  flagMothers(p, keepOrDrop); break;
216  case SelectCode::kNone:
217  ;
218  };
219  }
220  }
221  }
222  indices_.clear();
223  int counter = 0;
224  for(size_t i = 0; i < n; ++i) {
225  if(flags_[i] == keep) {
226  indices_.push_back(i);
227  flags_[i] = counter++;
228  } else
229  {
230  flags_[i]=-1; //set to invalid ref
231  }
232  }
233 
234  auto_ptr<GenParticleCollection> out(new GenParticleCollection);
236  out->reserve(counter);
237 
238  for(vector<size_t>::const_iterator i = indices_.begin(); i != indices_.end(); ++i) {
239  size_t index = *i;
240  const GenParticle & gen = (*src)[index];
241  const LeafCandidate & part = gen;
242  out->push_back(GenParticle(part));
243  GenParticle & newGen = out->back();
244  // The "daIndxs" and "moIndxs" keep a list of the keys for the mother/daughter
245  // parentage/descendency. In some cases, a circular referencing is encountered,
246  // which would result in an infinite loop. The list is checked to
247  // avoid this.
248  vector<size_t> daIndxs;
249  addDaughterRefs(daIndxs, newGen, outRef, gen.daughterRefVector());
250  vector<size_t> moIndxs;
251  addMotherRefs(moIndxs, newGen, outRef, gen.motherRefVector());
252  }
253 
254 
256  std::auto_ptr<edm::Association<reco::GenParticleCollection> > orig2new(new edm::Association<reco::GenParticleCollection>(oh ));
258  orig2newFiller.insert(src, flags_.begin(), flags_.end());
259  orig2newFiller.fill();
260  evt.put(orig2new);
261 
262 
263 }
264 
265 
266 void GenParticlePruner::addDaughterRefs(vector<size_t> & daIndxs,
267  GenParticle& newGen, GenParticleRefProd outRef,
268  const GenParticleRefVector& daughters) const {
269  for(GenParticleRefVector::const_iterator j = daughters.begin();
270  j != daughters.end(); ++j) {
271  GenParticleRef dau = *j;
272  if ( find(daIndxs.begin(), daIndxs.end(), dau.key()) == daIndxs.end() ) {
273  int idx = flags_[dau.key()];
274  daIndxs.push_back( dau.key() );
275  if(idx > 0) {
276  GenParticleRef newDau(outRef, static_cast<size_t>(idx));
277  newGen.addDaughter(newDau);
278  } else {
279  const GenParticleRefVector daus = dau->daughterRefVector();
280  if(daus.size()>0) {
281  addDaughterRefs(daIndxs, newGen, outRef, daus);
282  }
283  }
284  }
285  }
286 }
287 
288 
289 
290 void GenParticlePruner::addMotherRefs(vector<size_t> & moIndxs,
291  GenParticle& newGen, GenParticleRefProd outRef,
292  const GenParticleRefVector& mothers) const {
294  j != mothers.end(); ++j) {
295  GenParticleRef mom = *j;
296  if ( find(moIndxs.begin(), moIndxs.end(), mom.key()) == moIndxs.end() ) {
297  int idx = flags_[mom.key()];
298  moIndxs.push_back( mom.key() );
299  if(idx >= 0) {
300  GenParticleRef newMom(outRef, static_cast<size_t>(idx));
301  newGen.addMother(newMom);
302  } else {
303  const GenParticleRefVector moms = mom->motherRefVector();
304  if(moms.size()>0)
305  addMotherRefs(moIndxs, newGen, outRef, moms);
306  }
307  }
308  }
309 }
310 
312 
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_
tuple cfg
Definition: looper.py:237
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:446
#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 insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
void recursiveFlagMothers(size_t, const reco::GenParticleCollection &, int, std::vector< size_t > &)
const daughters & daughterRefVector() const
references to daughtes
GenParticlePruner(const edm::ParameterSet &)
key_type key() const
Accessor for product key.
Definition: Ref.h:266
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:113
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:133
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
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