CMS 3D CMS Logo

GenParticlePruner.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/EDProducer.h"
00002 #include "FWCore/ParameterSet/interface/InputTag.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00007 #include "PhysicsTools/Utilities/interface/StringCutObjectSelector.h"
00008 #include "FWCore/Utilities/interface/EDMException.h"
00009 #include "PhysicsTools/HepMCCandAlgos/interface/PdgEntryReplacer.h"
00010 
00011 namespace helper {
00012   struct SelectCode {
00013     enum KeepOrDrop { kKeep, kDrop };
00014     enum FlagDepth { kNone, kFirst, kAll };
00015     KeepOrDrop keepOrDrop_;
00016     FlagDepth daughtersDepth_, mothersDepth_;
00017     bool all_;
00018   };
00019 }
00020 
00021 class GenParticlePruner : public edm::EDProducer {
00022 public:
00023   GenParticlePruner(const edm::ParameterSet&);
00024 private:
00025   void produce(edm::Event&, const edm::EventSetup&);
00026   void beginJob(const edm::EventSetup&);
00027   edm::InputTag src_;
00028   int keepOrDropAll_;
00029   std::vector<std::string> selection_;
00030   std::vector<std::pair<StringCutObjectSelector<reco::GenParticle>, helper::SelectCode> > select_;
00031   std::vector<int> flags_;
00032   std::vector<size_t> indices_;
00033   void parse(const std::string & selection, helper::SelectCode & code, std::string & cut) const;
00034   void flagDaughters(const reco::GenParticle &, int); 
00035   void flagMothers(const reco::GenParticle &, int); 
00036   void recursiveFlagDaughters(size_t, const reco::GenParticleCollection &, int); 
00037   void recursiveFlagMothers(size_t, const reco::GenParticleCollection &, int); 
00038   void addDaughterRefs(reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
00039   void addMotherRefs(reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
00040 };
00041 
00042 using namespace edm;
00043 using namespace std;
00044 using namespace reco;
00045 
00046 const int keep = 1, drop = -1;
00047 
00048 void GenParticlePruner::parse(const std::string & selection, ::helper::SelectCode & code, std::string & cut) const {
00049   using namespace ::helper;
00050   size_t f =  selection.find_first_not_of(' ');
00051   size_t n = selection.size();
00052   string command;
00053   char c;
00054   for(; (c = selection[f]) != ' ' && f < n; ++f) {
00055     command.push_back(c);
00056   }
00057   if(command[0] == '+') {
00058     command.erase(0, 1);
00059     if(command[0] == '+') {
00060       command.erase(0, 1);
00061       code.mothersDepth_ = SelectCode::kAll;
00062     } else {
00063       code.mothersDepth_ = SelectCode::kFirst;
00064     }
00065   } else 
00066     code.mothersDepth_ = SelectCode::kNone;
00067 
00068   if(command[command.size() - 1] == '+') {
00069     command.erase(command.size() - 1);
00070     if(command[command.size()-1] == '+') {
00071       command.erase(command.size() - 1);
00072       code.daughtersDepth_ = SelectCode::kAll;
00073     } else {
00074       code.daughtersDepth_ = SelectCode::kFirst;
00075     }
00076   } else
00077     code.daughtersDepth_ = SelectCode::kNone;
00078 
00079   if(command == "keep") code.keepOrDrop_ = SelectCode::kKeep;
00080   else if(command == "drop") code.keepOrDrop_ = SelectCode::kDrop;
00081   else {
00082     throw Exception(errors::Configuration)
00083       << "invalid selection command: " << command << "\n" << endl;
00084   }
00085   for(; f < n; ++f) {
00086     if(selection[f] != ' ') break;
00087   }
00088   cut = string(selection, f);
00089   if(cut[0] == '*')
00090     cut = string(cut, 0, cut.find_first_of(' '));
00091   code.all_ = cut == "*";
00092 }
00093 
00094 GenParticlePruner::GenParticlePruner(const ParameterSet& cfg) :
00095   src_(cfg.getParameter<InputTag>("src")), keepOrDropAll_(drop),
00096   selection_(cfg.getParameter<vector<string> >("select")) {
00097   using namespace ::helper;
00098   produces<GenParticleCollection>();
00099 }
00100 
00101 void GenParticlePruner::beginJob(const EventSetup& es) {
00102   PdgEntryReplacer rep(es);
00103   for(vector<string>::const_iterator i = selection_.begin(); i != selection_.end(); ++i) {
00104     string cut;
00105     ::helper::SelectCode code;
00106     parse(*i, code, cut);
00107     if(code.all_) {
00108       if(i != selection_.begin()) 
00109         throw Exception(errors::Configuration)
00110           << "selections \"keep *\" and \"drop *\" can be used only as first options. Here used in position # " 
00111           << (i - selection_.begin()) + 1 << "\n" << endl;
00112       switch(code.keepOrDrop_) {
00113           case ::helper::SelectCode::kDrop :
00114               keepOrDropAll_ = drop; break;
00115           case ::helper::SelectCode::kKeep :
00116               keepOrDropAll_ = keep; 
00117       };
00118     } else {
00119       cut = rep.replace(cut);
00120       select_.push_back(make_pair(StringCutObjectSelector<GenParticle>(cut), code));
00121     }
00122   }
00123 }
00124 
00125 void GenParticlePruner::flagDaughters(const reco::GenParticle & gen, int keepOrDrop) {
00126   GenParticleRefVector daughters = gen.daughterRefVector();
00127   for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i) 
00128     flags_[i->key()] = keepOrDrop;
00129 }
00130 
00131 void GenParticlePruner::flagMothers(const reco::GenParticle & gen, int keepOrDrop) {
00132   GenParticleRefVector mothers = gen.motherRefVector();
00133   for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i) 
00134     flags_[i->key()] = keepOrDrop;
00135 }
00136 
00137 void GenParticlePruner::recursiveFlagDaughters(size_t index, const reco::GenParticleCollection & src, int keepOrDrop) {
00138   GenParticleRefVector daughters = src[index].daughterRefVector();
00139   for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i) {
00140     index = i->key();
00141     flags_[index] = keepOrDrop;
00142     recursiveFlagDaughters(index, src, keepOrDrop);
00143   }
00144 }
00145 
00146 void GenParticlePruner::recursiveFlagMothers(size_t index, const reco::GenParticleCollection & src, int keepOrDrop) {
00147   GenParticleRefVector mothers = src[index].motherRefVector();
00148   for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i) {
00149     index = i->key();
00150     flags_[index] = keepOrDrop;
00151     recursiveFlagMothers(index, src, keepOrDrop);
00152   }
00153 }
00154 
00155 void GenParticlePruner::produce(Event& evt, const EventSetup&) {
00156   using namespace ::helper;
00157   Handle<GenParticleCollection> src;
00158   evt.getByLabel(src_, src);
00159   const size_t n = src->size();
00160   flags_.clear();
00161   flags_.resize(n, keepOrDropAll_);
00162   for(size_t j = 0; j < select_.size(); ++j) { 
00163     const pair<StringCutObjectSelector<GenParticle>, SelectCode> & sel = select_[j];
00164     SelectCode code = sel.second;
00165     const StringCutObjectSelector<GenParticle> & cut = sel.first;
00166     for(size_t i = 0; i < n; ++i) {
00167       const GenParticle & p = (*src)[i];
00168       if(cut(p)) {
00169         int keepOrDrop = keep;
00170         switch(code.keepOrDrop_) {
00171         case SelectCode::kKeep:
00172           keepOrDrop = keep; break;
00173         case SelectCode::kDrop:
00174           keepOrDrop = drop; 
00175         };
00176         flags_[i] = keepOrDrop;
00177         switch(code.daughtersDepth_) {
00178         case SelectCode::kAll : 
00179           recursiveFlagDaughters(i, *src, keepOrDrop); break;
00180         case SelectCode::kFirst :
00181           flagDaughters(p, keepOrDrop); break;
00182         case SelectCode::kNone:
00183           ;
00184         };
00185         switch(code.mothersDepth_) {
00186         case SelectCode::kAll :
00187           recursiveFlagMothers(i, *src, keepOrDrop); break;
00188         case SelectCode::kFirst :
00189           flagMothers(p, keepOrDrop); break;
00190         case SelectCode::kNone:
00191           ;
00192         };
00193       }
00194     }
00195   }
00196   indices_.clear();
00197   int counter = 0;
00198   for(size_t i = 0; i < n; ++i) {
00199     if(flags_[i] == keep) {
00200       indices_.push_back(i);
00201       flags_[i] = counter++;
00202     }
00203   }
00204 
00205   auto_ptr<GenParticleCollection> out(new GenParticleCollection);
00206   GenParticleRefProd outRef = evt.getRefBeforePut<GenParticleCollection>();
00207   out->reserve(counter);
00208   for(vector<size_t>::const_iterator i = indices_.begin(); i != indices_.end(); ++i) {
00209     size_t index = *i;
00210     const GenParticle & gen = (*src)[index];
00211     const Particle & part = gen;
00212     out->push_back(GenParticle(part));
00213     GenParticle & newGen = out->back();
00214     addDaughterRefs(newGen, outRef, gen.daughterRefVector());
00215     addMotherRefs(newGen, outRef, gen.motherRefVector());
00216   }
00217 
00218   evt.put(out);
00219 }
00220 
00221 
00222 void GenParticlePruner::addDaughterRefs(GenParticle& newGen, GenParticleRefProd outRef, 
00223                                         const GenParticleRefVector& daughters) const {
00224   for(GenParticleRefVector::const_iterator j = daughters.begin();
00225       j != daughters.end(); ++j) {
00226     GenParticleRef dau = *j;
00227     int idx = flags_[dau.key()];
00228     if(idx > 0) {
00229       GenParticleRef newDau(outRef, static_cast<size_t>(idx));
00230       newGen.addDaughter(newDau);
00231     } else {
00232       const GenParticleRefVector daus = dau->daughterRefVector();
00233       if(daus.size()>0)
00234         addDaughterRefs(newGen, outRef, daus);
00235     }
00236   }
00237 }
00238 
00239 void GenParticlePruner::addMotherRefs(GenParticle& newGen, GenParticleRefProd outRef, 
00240                                       const GenParticleRefVector& mothers) const {
00241   for(GenParticleRefVector::const_iterator j = mothers.begin();
00242       j != mothers.end(); ++j) {
00243     GenParticleRef mom = *j;
00244     int idx = flags_[mom.key()];
00245     if(idx > 0) {
00246       GenParticleRef newMom(outRef, static_cast<size_t>(idx));
00247       newGen.addMother(newMom);
00248     } else {
00249       const GenParticleRefVector moms = mom->motherRefVector();
00250       if(moms.size()>0)
00251         addMotherRefs(newGen, outRef, moms);
00252     }
00253   }
00254 }
00255 
00256 #include "FWCore/Framework/interface/MakerMacros.h"
00257 
00258 DEFINE_FWK_MODULE(GenParticlePruner);

Generated on Tue Jun 9 17:41:08 2009 for CMSSW by  doxygen 1.5.4