CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/PhysicsTools/HepMCCandAlgos/plugins/GenPlusSimParticleProducer.cc

Go to the documentation of this file.
00001 
00022 #include "FWCore/Framework/interface/Frameworkfwd.h"
00023 #include "FWCore/Framework/interface/EDProducer.h"
00024 #include "FWCore/Framework/interface/Event.h"
00025 #include "FWCore/Framework/interface/MakerMacros.h"
00026 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00027 
00028 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00029 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00030 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00031 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00032 
00033 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00034 #include "SimGeneral/HepPDTRecord/interface/PdtEntry.h"
00035 
00036 #include <ext/algorithm>
00037 
00038 namespace pat {
00039 class GenPlusSimParticleProducer : public edm::EDProducer {
00040 public:
00041   explicit GenPlusSimParticleProducer(const edm::ParameterSet&);
00042   ~GenPlusSimParticleProducer() {}
00043 
00044 private:
00045   virtual void produce(edm::Event&, const edm::EventSetup&);
00046   virtual void endJob() {}
00047 
00048   bool firstEvent_;
00049   edm::InputTag src_;
00050   int setStatus_;
00051   std::set<int>         pdgIds_; // these are the ones we really use
00052   std::vector<PdtEntry> pdts_;   // these are needed before we get the EventSetup
00053 
00054   typedef StringCutObjectSelector<reco::GenParticle> StrFilter;
00055   std::auto_ptr<StrFilter> filter_;
00056 
00058   edm::InputTag genParticles_;
00059 
00061 
00074   void addGenParticle(const SimTrack &stMom,
00075                       const SimTrack &stDau,
00076                       unsigned int momidx,
00077                       const edm::SimTrackContainer &simtks,
00078                       const edm::SimVertexContainer &simvtxs,
00079                       reco::GenParticleCollection &mergedGens,
00080                       const reco::GenParticleRefProd ref,
00081                       std::vector<int> &genBarcodes,
00082                       bool &barcodesAreSorted ) const ;
00083   struct LessById {
00084     bool operator()(const SimTrack &tk1, const SimTrack &tk2) const { return tk1.trackId() < tk2.trackId(); }
00085     bool operator()(const SimTrack &tk1, unsigned int    id ) const { return tk1.trackId() < id;            }
00086     bool operator()(unsigned int     id, const SimTrack &tk2) const { return id            < tk2.trackId(); }
00087   };
00088   
00089 };
00090 }
00091 
00092 using namespace std;
00093 using namespace edm;
00094 using namespace reco;
00095 using pat::GenPlusSimParticleProducer;
00096 
00097 GenPlusSimParticleProducer::GenPlusSimParticleProducer(const ParameterSet& cfg) :
00098   firstEvent_(true),
00099   src_(cfg.getParameter<InputTag>("src")),            // source sim tracks & vertices
00100   setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
00101   genParticles_(cfg.getParameter<InputTag>("genParticles")) // get the genParticles to add GEANT particles to
00102 {
00103   // Possibly allow a list of particle types
00104   if (cfg.exists("particleTypes")) {
00105     pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
00106   }
00107   
00108   // Possibly allow a string cut
00109   if (cfg.existsAs<string>("filter")) {
00110     string filter = cfg.getParameter<string>("filter");
00111     if (!filter.empty()) {
00112       filter_ = auto_ptr<StrFilter>(new StrFilter(filter));
00113     }
00114   }
00115   produces<GenParticleCollection>();
00116   produces<vector<int> >();
00117 }
00118 
00119 void GenPlusSimParticleProducer::addGenParticle( const SimTrack &stMom,
00120                                                     const SimTrack &stDau,
00121                                                     unsigned int momidx,
00122                                                     const SimTrackContainer &simtracksSorted, 
00123                                                     const SimVertexContainer &simvertices, 
00124                                                     reco::GenParticleCollection &mergedGens,
00125                                                     const GenParticleRefProd ref,
00126                                                     std::vector<int> &genBarcodes,
00127                                                     bool &barcodesAreSorted) const 
00128 {
00129   // Make the genParticle for stDau and add it to the new collection and update the parent-child relationship
00130   // Make up a GenParticleCandidate from the GEANT track info.
00131   int charge = static_cast<int>(stDau.charge());
00132   Particle::LorentzVector p4 = stDau.momentum();
00133   Particle::Point vtx; // = (0,0,0) by default
00134   if (!stDau.noVertex()) vtx = simvertices[stDau.vertIndex()].position();
00135   GenParticle genp(charge, p4, vtx, stDau.type(), setStatus_, true);
00136   
00137   // Maybe apply filter on the particle
00138   if (filter_.get() != 0) {
00139     if (!(*filter_)(genp)) return;
00140   }
00141   
00142   reco::GenParticleRef parentRef(ref, momidx);
00143   genp.addMother(parentRef);
00144   mergedGens.push_back(genp);
00145   // get the index for the daughter just added
00146   unsigned int dauidx = mergedGens.size()-1;
00147 
00148   // update add daughter relationship
00149   reco::GenParticle & cand = mergedGens[ momidx ];
00150   cand.addDaughter( GenParticleRef( ref, dauidx) );
00151 
00152   //look for simtrack daughters of stDau to see if we need to recur further down the chain
00153   
00154   for (SimTrackContainer::const_iterator isimtrk = simtracksSorted.begin();
00155        isimtrk != simtracksSorted.end(); ++isimtrk) {
00156     if (!isimtrk->noVertex()) {
00157       // Pick the vertex (isimtrk.vertIndex() is really an index)
00158       const SimVertex &vtx = simvertices[isimtrk->vertIndex()];
00159       
00160       // Check if the vertex has a parent track (otherwise, we're lost)
00161       if (!vtx.noParent()) {
00162         
00163         // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it 
00164         unsigned int idx = vtx.parentIndex();
00165         SimTrackContainer::const_iterator it = std::lower_bound(simtracksSorted.begin(), simtracksSorted.end(), idx, LessById());
00166         if ((it != simtracksSorted.end()) && (it->trackId() == idx)) {
00167           if (it->trackId()==stDau.trackId()) {
00168             //need the genparticle index of stDau which is dauidx
00169             addGenParticle(stDau, *isimtrk, dauidx, simtracksSorted, simvertices, mergedGens, ref, genBarcodes, barcodesAreSorted);
00170           }
00171         }
00172       }
00173     }
00174   }
00175 }
00176 
00177 void GenPlusSimParticleProducer::produce(Event& event,
00178                                             const EventSetup& iSetup) {
00179   if (firstEvent_){
00180     if (!pdts_.empty()) {
00181       pdgIds_.clear();
00182       for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
00183         itp->setup(iSetup); // decode string->pdgId and vice-versa                                                                                                      
00184         pdgIds_.insert(std::abs(itp->pdgId()));
00185       }
00186       pdts_.clear();
00187     }
00188     firstEvent_ = false;
00189   }
00190 
00191   // Simulated tracks (i.e. GEANT particles).
00192   Handle<SimTrackContainer> simtracks;
00193   event.getByLabel(src_, simtracks);
00194   
00195   // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
00196   std::auto_ptr<SimTrackContainer> simtracksTmp;
00197   const SimTrackContainer * simtracksSorted = &* simtracks;
00198   if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
00199     simtracksTmp.reset(new SimTrackContainer(*simtracks));
00200     std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
00201     simtracksSorted = &* simtracksTmp;
00202   }
00203   
00204   // Get the associated vertices
00205   Handle<SimVertexContainer> simvertices;
00206   event.getByLabel(src_, simvertices);
00207   
00208   // Get the GenParticles and barcodes, if needed to set mother and daughter links
00209   Handle<GenParticleCollection> gens;
00210   Handle<std::vector<int> > genBarcodes;
00211   bool barcodesAreSorted = true;
00212   event.getByLabel(genParticles_, gens);
00213   event.getByLabel(genParticles_, genBarcodes);
00214   if (gens->size() != genBarcodes->size()) throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
00215   
00216   // make the output collection
00217   auto_ptr<GenParticleCollection> candsPtr(new GenParticleCollection);
00218   GenParticleCollection & cands = * candsPtr;
00219   
00220   const GenParticleRefProd ref = event.getRefBeforePut<GenParticleCollection>();
00221   
00222   // add the original genParticles to the merged output list
00223   for( size_t i = 0; i < gens->size(); ++ i ) {
00224     reco::GenParticle cand((*gens)[i]);
00225     cands.push_back(cand);
00226   }
00227   
00228   // make new barcodes vector and fill it with the original list
00229   auto_ptr<vector<int> > newGenBarcodes( new vector<int>() );
00230   for (unsigned int i = 0; i < genBarcodes->size(); ++i) {
00231     newGenBarcodes->push_back((*genBarcodes)[i]);
00232   }
00233   barcodesAreSorted = __gnu_cxx::is_sorted(newGenBarcodes->begin(), newGenBarcodes->end());
00234   
00235   for( size_t i = 0; i < cands.size(); ++ i ) {
00236     reco::GenParticle & cand = cands[ i ];
00237     size_t nDaus = cand.numberOfDaughters();
00238     GenParticleRefVector daus = cand.daughterRefVector();
00239     cand.resetDaughters( ref.id() );
00240     for ( size_t d = 0; d < nDaus; ++d) {
00241       cand.addDaughter( GenParticleRef( ref, daus[d].key() ) );
00242     }
00243     
00244     size_t nMoms = cand.numberOfMothers();
00245     GenParticleRefVector moms = cand.motherRefVector();
00246     cand.resetMothers( ref.id() );
00247     for ( size_t m = 0; m < nMoms; ++m) {
00248       cand.addMother( GenParticleRef( ref, moms[m].key() ) );
00249     }
00250   }
00251   
00252   for (SimTrackContainer::const_iterator isimtrk = simtracks->begin();
00253        isimtrk != simtracks->end(); ++isimtrk) {
00254     
00255     // Skip PYTHIA tracks.
00256     if (isimtrk->genpartIndex() != -1) continue; 
00257     
00258     // Maybe apply the PdgId filter
00259     if (!pdgIds_.empty()) { // if we have a filter on pdg ids
00260       if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end()) continue;
00261     }
00262     
00263     // find simtrack that has a genParticle match to its parent
00264     // Look at the production vertex. If there is no vertex, I can do nothing...
00265     if (!isimtrk->noVertex()) {
00266       
00267       // Pick the vertex (isimtrk.vertIndex() is really an index)
00268       const SimVertex &vtx = (*simvertices)[isimtrk->vertIndex()];
00269       
00270       // Check if the vertex has a parent track (otherwise, we're lost)
00271       if (!vtx.noParent()) {
00272         
00273         // Now note that vtx.parentIndex() is NOT an index, it's a track id, so I have to search for it 
00274         unsigned int idx = vtx.parentIndex();
00275         SimTrackContainer::const_iterator it = std::lower_bound(simtracksSorted->begin(), simtracksSorted->end(), idx, LessById());
00276         if ((it != simtracksSorted->end()) && (it->trackId() == idx)) { //it is the parent sim track
00277           if (it->genpartIndex() != -1) {
00278             std::vector<int>::const_iterator itIndex;
00279             if (barcodesAreSorted) {
00280               itIndex = std::lower_bound(genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
00281             } else {
00282               itIndex = std::find(       genBarcodes->begin(), genBarcodes->end(), it->genpartIndex());
00283             }
00284             
00285             // Check that I found something
00286             // I need to check '*itIndex == it->genpartIndex()' because lower_bound just finds the right spot for an item in a sorted list, not the item
00287             if ((itIndex != genBarcodes->end()) && (*itIndex == it->genpartIndex())) {
00288               // Ok, I'll make the genParticle for st and add it to the new collection updating the map and parent-child relationship
00289               // pass the mother and daughter sim tracks and the mother genParticle to method to create the daughter genParticle and recur
00290               unsigned int momidx = itIndex - genBarcodes->begin();
00291               addGenParticle(*it, *isimtrk, momidx, *simtracksSorted, *simvertices, *candsPtr, ref, *newGenBarcodes, barcodesAreSorted);
00292             }
00293           }
00294         }
00295       }
00296     }   
00297   }
00298   
00299   event.put(candsPtr);
00300   event.put(newGenBarcodes);
00301 }
00302 
00303 DEFINE_FWK_MODULE(GenPlusSimParticleProducer);