CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/PhysicsTools/HepMCCandAlgos/plugins/GenParticlePruner.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/EDProducer.h"
00002 #include "FWCore/Utilities/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 "CommonTools/Utils/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   bool firstEvent_;
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, std::vector<size_t> &); 
00037   void recursiveFlagMothers(size_t, const reco::GenParticleCollection &, int, std::vector<size_t> &); 
00038   void addDaughterRefs(std::vector<size_t> &, reco::GenParticle&, reco::GenParticleRefProd, const reco::GenParticleRefVector&) const;
00039   void addMotherRefs(std::vector<size_t> &, 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   firstEvent_(true),
00096   src_(cfg.getParameter<InputTag>("src")), keepOrDropAll_(drop),
00097   selection_(cfg.getParameter<vector<string> >("select")) {
00098   using namespace ::helper;
00099   produces<GenParticleCollection>();
00100 }
00101 
00102 void GenParticlePruner::flagDaughters(const reco::GenParticle & gen, int keepOrDrop) {
00103   GenParticleRefVector daughters = gen.daughterRefVector();
00104   for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i) 
00105     flags_[i->key()] = keepOrDrop;
00106 }
00107 
00108 void GenParticlePruner::flagMothers(const reco::GenParticle & gen, int keepOrDrop) {
00109   GenParticleRefVector mothers = gen.motherRefVector();
00110   for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i) 
00111     flags_[i->key()] = keepOrDrop;
00112 }
00113 
00114 void GenParticlePruner::recursiveFlagDaughters(size_t index, const reco::GenParticleCollection & src, int keepOrDrop,
00115                                                std::vector<size_t> & allIndices ) {
00116   GenParticleRefVector daughters = src[index].daughterRefVector();
00117   // avoid infinite recursion if the daughters are set to "this" particle.
00118   size_t cachedIndex = index;
00119   for(GenParticleRefVector::const_iterator i = daughters.begin(); i != daughters.end(); ++i) {
00120     index = i->key();
00121     // To also avoid infinite recursion if a "loop" is found in the daughter list,
00122     // check to make sure the index hasn't already been added. 
00123     if ( find( allIndices.begin(), allIndices.end(), index ) == allIndices.end() ) {
00124       allIndices.push_back( index );
00125       if ( cachedIndex != index ) {
00126         flags_[index] = keepOrDrop;
00127         recursiveFlagDaughters(index, src, keepOrDrop, allIndices);
00128       }
00129     }
00130   }
00131 }
00132 
00133 void GenParticlePruner::recursiveFlagMothers(size_t index, const reco::GenParticleCollection & src, int keepOrDrop,
00134                                              std::vector<size_t> & allIndices ) {
00135   GenParticleRefVector mothers = src[index].motherRefVector();
00136   // avoid infinite recursion if the mothers are set to "this" particle.
00137   size_t cachedIndex = index;
00138   for(GenParticleRefVector::const_iterator i = mothers.begin(); i != mothers.end(); ++i) {
00139     index = i->key();
00140     // To also avoid infinite recursion if a "loop" is found in the daughter list,
00141     // check to make sure the index hasn't already been added. 
00142     if ( find( allIndices.begin(), allIndices.end(), index ) == allIndices.end() ) {
00143       allIndices.push_back( index );
00144       if ( cachedIndex != index ) {
00145         flags_[index] = keepOrDrop;
00146         recursiveFlagMothers(index, src, keepOrDrop, allIndices);
00147       }
00148     }
00149   }
00150 }
00151 
00152 void GenParticlePruner::produce(Event& evt, const EventSetup& es) {
00153   if (firstEvent_) {
00154     PdgEntryReplacer rep(es);
00155     for(vector<string>::const_iterator i = selection_.begin(); i != selection_.end(); ++i) {
00156       string cut;
00157       ::helper::SelectCode code;
00158       parse(*i, code, cut);
00159       if(code.all_) {
00160         if(i != selection_.begin())
00161           throw Exception(errors::Configuration)
00162             << "selections \"keep *\" and \"drop *\" can be used only as first options. Here used in position # "
00163             << (i - selection_.begin()) + 1 << "\n" << endl;
00164         switch(code.keepOrDrop_) {
00165         case ::helper::SelectCode::kDrop :
00166           keepOrDropAll_ = drop; break;
00167         case ::helper::SelectCode::kKeep :
00168           keepOrDropAll_ = keep;
00169         };
00170       } else {
00171         cut = rep.replace(cut);
00172         select_.push_back(make_pair(StringCutObjectSelector<GenParticle>(cut), code));
00173       }
00174     }
00175     firstEvent_ = false;
00176   }
00177 
00178   using namespace ::helper;
00179   Handle<GenParticleCollection> src;
00180   evt.getByLabel(src_, src);
00181   const size_t n = src->size();
00182   flags_.clear();
00183   flags_.resize(n, keepOrDropAll_);
00184   for(size_t j = 0; j < select_.size(); ++j) { 
00185     const pair<StringCutObjectSelector<GenParticle>, SelectCode> & sel = select_[j];
00186     SelectCode code = sel.second;
00187     const StringCutObjectSelector<GenParticle> & cut = sel.first;
00188     for(size_t i = 0; i < n; ++i) {
00189       const GenParticle & p = (*src)[i];
00190       if(cut(p)) {
00191         int keepOrDrop = keep;
00192         switch(code.keepOrDrop_) {
00193         case SelectCode::kKeep:
00194           keepOrDrop = keep; break;
00195         case SelectCode::kDrop:
00196           keepOrDrop = drop; 
00197         };
00198         flags_[i] = keepOrDrop;
00199         std::vector<size_t> allIndicesDa;
00200         std::vector<size_t> allIndicesMo;
00201         switch(code.daughtersDepth_) {
00202         case SelectCode::kAll : 
00203           recursiveFlagDaughters(i, *src, keepOrDrop, allIndicesDa); break;
00204         case SelectCode::kFirst :
00205           flagDaughters(p, keepOrDrop); break;
00206         case SelectCode::kNone:
00207           ;
00208         };
00209         switch(code.mothersDepth_) {
00210         case SelectCode::kAll :
00211           recursiveFlagMothers(i, *src, keepOrDrop, allIndicesMo); break;
00212         case SelectCode::kFirst :
00213           flagMothers(p, keepOrDrop); break;
00214         case SelectCode::kNone:
00215           ;
00216         };
00217       }
00218     }
00219   }
00220   indices_.clear();
00221   int counter = 0;
00222   for(size_t i = 0; i < n; ++i) {
00223     if(flags_[i] == keep) {
00224       indices_.push_back(i);
00225       flags_[i] = counter++;
00226     }
00227   }
00228 
00229   auto_ptr<GenParticleCollection> out(new GenParticleCollection);
00230   GenParticleRefProd outRef = evt.getRefBeforePut<GenParticleCollection>();
00231   out->reserve(counter);
00232   for(vector<size_t>::const_iterator i = indices_.begin(); i != indices_.end(); ++i) {
00233     size_t index = *i;
00234     const GenParticle & gen = (*src)[index];
00235     const LeafCandidate & part = gen;
00236     out->push_back(GenParticle(part));
00237     GenParticle & newGen = out->back();
00238     // The "daIndxs" and "moIndxs" keep a list of the keys for the mother/daughter
00239     // parentage/descendency. In some cases, a circular referencing is encountered,
00240     // which would result in an infinite loop. The list is checked to
00241     // avoid this. 
00242     vector<size_t> daIndxs;
00243     addDaughterRefs(daIndxs, newGen, outRef, gen.daughterRefVector());
00244     vector<size_t> moIndxs;
00245     addMotherRefs(moIndxs, newGen, outRef, gen.motherRefVector());
00246   }
00247 
00248   evt.put(out);
00249 }
00250 
00251 
00252 void GenParticlePruner::addDaughterRefs(vector<size_t> & daIndxs, 
00253                                         GenParticle& newGen, GenParticleRefProd outRef, 
00254                                         const GenParticleRefVector& daughters) const {
00255   for(GenParticleRefVector::const_iterator j = daughters.begin();
00256       j != daughters.end(); ++j) {
00257     GenParticleRef dau = *j;
00258     if ( find(daIndxs.begin(), daIndxs.end(), dau.key()) == daIndxs.end() ) {
00259       int idx = flags_[dau.key()];
00260       daIndxs.push_back( dau.key() );
00261       if(idx > 0) {
00262         GenParticleRef newDau(outRef, static_cast<size_t>(idx));
00263         newGen.addDaughter(newDau);
00264       } else {
00265         const GenParticleRefVector daus = dau->daughterRefVector();
00266         if(daus.size()>0) {
00267           addDaughterRefs(daIndxs, newGen, outRef, daus);
00268         }
00269       }
00270     }
00271   }
00272 }
00273 
00274 
00275 
00276 void GenParticlePruner::addMotherRefs(vector<size_t> & moIndxs,
00277                                       GenParticle& newGen, GenParticleRefProd outRef, 
00278                                       const GenParticleRefVector& mothers) const {
00279   for(GenParticleRefVector::const_iterator j = mothers.begin();
00280       j != mothers.end(); ++j) {
00281     GenParticleRef mom = *j;
00282     if ( find(moIndxs.begin(), moIndxs.end(), mom.key()) == moIndxs.end() ) {
00283       int idx = flags_[mom.key()];
00284       moIndxs.push_back( mom.key() );
00285       if(idx > 0) {
00286         GenParticleRef newMom(outRef, static_cast<size_t>(idx));
00287         newGen.addMother(newMom);
00288       } else {
00289         const GenParticleRefVector moms = mom->motherRefVector();
00290         if(moms.size()>0)
00291           addMotherRefs(moIndxs, newGen, outRef, moms);
00292       }
00293     }
00294   }
00295 }
00296 
00297 #include "FWCore/Framework/interface/MakerMacros.h"
00298 
00299 DEFINE_FWK_MODULE(GenParticlePruner);