CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/PhysicsTools/PatAlgos/plugins/PATGenCandsFromSimTracksProducer.cc

Go to the documentation of this file.
00001 
00012 #include "FWCore/Framework/interface/Frameworkfwd.h"
00013 #include "FWCore/Framework/interface/EDProducer.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/MakerMacros.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 
00018 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00019 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00020 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00021 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00022 
00023 #include "CommonTools/Utils/interface/StringCutObjectSelector.h"
00024 #include "SimGeneral/HepPDTRecord/interface/PdtEntry.h"
00025 
00026 #include <ext/algorithm>
00027 
00028 namespace pat {
00029 class PATGenCandsFromSimTracksProducer : public edm::EDProducer {
00030 public:
00031   explicit PATGenCandsFromSimTracksProducer(const edm::ParameterSet&);
00032   ~PATGenCandsFromSimTracksProducer() {}
00033 
00034 private:
00035   virtual void produce(edm::Event&, const edm::EventSetup&);
00036   virtual void endJob() {}
00037 
00038   bool firstEvent_;
00039   edm::InputTag src_;
00040   int setStatus_;
00041   std::set<int>         pdgIds_; // these are the ones we really use
00042   std::vector<PdtEntry> pdts_;   // these are needed before we get the EventSetup
00043   std::set<int>         motherPdgIds_; // these are the ones we really use
00044   std::vector<PdtEntry> motherPdts_;   // these are needed before we get the EventSetup
00045 
00046   typedef StringCutObjectSelector<reco::GenParticle> StrFilter;
00047   std::auto_ptr<StrFilter> filter_;
00048 
00050   bool makeMotherLink_;
00052   bool writeAncestors_;
00053 
00055   edm::InputTag genParticles_;
00056 
00058   struct GlobalContext {
00059       GlobalContext(const edm::SimTrackContainer &simtks1,
00060                     const edm::SimVertexContainer &simvtxs1,
00061                     const edm::Handle<reco::GenParticleCollection> &gens1,
00062                     const edm::Handle<std::vector<int> >           &genBarcodes1,
00063                     bool                                            barcodesAreSorted1,
00064                     reco::GenParticleCollection                     & output1,
00065                     const edm::RefProd<reco::GenParticleCollection> & refprod1) :
00066           simtks(simtks1), simvtxs(simvtxs1), 
00067           gens(gens1), genBarcodes(genBarcodes1), barcodesAreSorted(barcodesAreSorted1),
00068           output(output1), refprod(refprod1), simTksProcessed() {}
00069       // GEANT info
00070       const edm::SimTrackContainer &simtks;
00071       const edm::SimVertexContainer &simvtxs;
00072       // PYTHIA info
00073       const edm::Handle<reco::GenParticleCollection> &gens;
00074       const edm::Handle<std::vector<int> >           &genBarcodes;
00075       bool                                            barcodesAreSorted;
00076       // MY OUTPUT info
00077       reco::GenParticleCollection                     & output;
00078       const edm::RefProd<reco::GenParticleCollection> & refprod;
00079       // BOOK-KEEPING
00080       std::map<unsigned int,int> simTksProcessed; // key = sim track id; 
00081                                          // val = 0: not processed; 
00082                                          //       i>0:  (index+1) in my output
00083                                          //       i<0: -(index+1) in pythia [NOT USED]
00084   };
00085 
00087   const SimTrack * findGeantMother(const SimTrack &tk, const GlobalContext &g) const ;
00093   edm::Ref<reco::GenParticleCollection> findRef(const SimTrack &tk, GlobalContext &g) const ;
00094 
00096   edm::Ref<reco::GenParticleCollection> generatorRef_(const SimTrack &tk, const GlobalContext &g) const ;
00098   reco::GenParticle makeGenParticle_(const SimTrack &tk, const edm::Ref<reco::GenParticleCollection> & mother, const GlobalContext &g) const ;
00099   
00100 
00101 
00102   struct LessById {
00103     bool operator()(const SimTrack &tk1, const SimTrack &tk2) const { return tk1.trackId() < tk2.trackId(); }
00104     bool operator()(const SimTrack &tk1, unsigned int    id ) const { return tk1.trackId() < id;            }
00105     bool operator()(unsigned int     id, const SimTrack &tk2) const { return id            < tk2.trackId(); }
00106   };
00107    
00108 };
00109 }
00110 
00111 using namespace std;
00112 using namespace edm;
00113 using namespace reco;
00114 using pat::PATGenCandsFromSimTracksProducer;
00115 
00116 PATGenCandsFromSimTracksProducer::PATGenCandsFromSimTracksProducer(const ParameterSet& cfg) :
00117   firstEvent_(true),
00118   src_(cfg.getParameter<InputTag>("src")),            // source sim tracks & vertices
00119   setStatus_(cfg.getParameter<int32_t>("setStatus")), // set status of GenParticle to this code
00120   makeMotherLink_(cfg.existsAs<bool>("makeMotherLink") ? cfg.getParameter<bool>("makeMotherLink") : false),
00121   writeAncestors_(cfg.existsAs<bool>("writeAncestors") ? cfg.getParameter<bool>("writeAncestors") : false),
00122   genParticles_(makeMotherLink_ ? cfg.getParameter<InputTag>("genParticles") : edm::InputTag())
00123 {
00124     // Possibly allow a list of particle types
00125     if (cfg.exists("particleTypes")) {
00126         pdts_ = cfg.getParameter<vector<PdtEntry> >("particleTypes");
00127     }
00128     if (cfg.exists("motherTypes")) {
00129         motherPdts_ = cfg.getParameter<vector<PdtEntry> >("motherTypes");
00130     }
00131 
00132     // Possibly allow a string cut
00133     if (cfg.existsAs<string>("filter")) {
00134         string filter = cfg.getParameter<string>("filter");
00135         if (!filter.empty()) {
00136             filter_ = auto_ptr<StrFilter>(new StrFilter(filter));
00137         }
00138     }
00139 
00140     if (writeAncestors_ && !makeMotherLink_) {
00141         edm::LogWarning("Configuration") << "PATGenCandsFromSimTracksProducer: " << 
00142             "you have set 'writeAncestors' to 'true' and 'makeMotherLink' to false;" <<
00143             "GEANT particles with generator level (e.g.PYHIA) mothers won't have mother links.\n";
00144     }
00145     produces<GenParticleCollection>();
00146 }
00147 
00148 const SimTrack * 
00149 PATGenCandsFromSimTracksProducer::findGeantMother(const SimTrack &tk, const GlobalContext &g) const {
00150    assert(tk.genpartIndex() == -1); // MUST NOT be called with a PYTHIA track
00151    if (!tk.noVertex()) {
00152        const SimVertex &vtx = g.simvtxs[tk.vertIndex()];
00153        if (!vtx.noParent()) {
00154            unsigned int idx = vtx.parentIndex();
00155            SimTrackContainer::const_iterator it = std::lower_bound(g.simtks.begin(), g.simtks.end(), idx, LessById());
00156            if ((it != g.simtks.end()) && (it->trackId() == idx)) {
00157                 return &*it;
00158            }
00159        }
00160    }
00161    return 0;
00162 }
00163 
00164 edm::Ref<reco::GenParticleCollection> 
00165 PATGenCandsFromSimTracksProducer::findRef(const SimTrack &tk, GlobalContext &g) const {
00166     if (tk.genpartIndex() != -1) return makeMotherLink_ ? generatorRef_(tk, g) : edm::Ref<reco::GenParticleCollection>();
00167     const SimTrack * simMother = findGeantMother(tk, g);
00168 
00169     edm::Ref<reco::GenParticleCollection> motherRef;
00170     if (simMother != 0) motherRef = findRef(*simMother,g);
00171 
00172     if (writeAncestors_) {
00173         // If writing ancestors, I need to serialize myself, and then to return a ref to me
00174         // But first check if I've already been serialized
00175         std::map<unsigned int,int>::const_iterator it = g.simTksProcessed.find(tk.trackId());
00176         if (it != g.simTksProcessed.end()) {
00177             // just return a ref to it
00178             assert(it->second > 0);
00179             return edm::Ref<reco::GenParticleCollection>(g.refprod, (it->second) - 1);
00180         } else {
00181             // make genParticle, save, update the map, and return ref to myself
00182             reco::GenParticle p = makeGenParticle_(tk, motherRef, g); 
00183             g.output.push_back(p);
00184             g.simTksProcessed[tk.trackId()] = g.output.size();
00185             return edm::Ref<reco::GenParticleCollection>(g.refprod, g.output.size()-1 );
00186         }
00187     } else {
00188         // Otherwise, I just return a ref to my mum
00189         return motherRef; 
00190     }
00191 }
00192 
00193 edm::Ref<reco::GenParticleCollection> 
00194 PATGenCandsFromSimTracksProducer::generatorRef_(const SimTrack &st, const GlobalContext &g) const {
00195     assert(st.genpartIndex() != -1);        
00196     // Note that st.genpartIndex() is the barcode, not the index within GenParticleCollection, so I have to search the particle
00197     std::vector<int>::const_iterator it;
00198     if (g.barcodesAreSorted) {
00199         it = std::lower_bound(g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
00200     } else {
00201         it = std::find(       g.genBarcodes->begin(), g.genBarcodes->end(), st.genpartIndex());
00202     }
00203 
00204     // Check that I found something
00205     // I need to check '*it == st.genpartIndex()' because lower_bound just finds the right spot for an item in a sorted list, not the item
00206     if ((it != g.genBarcodes->end()) && (*it == st.genpartIndex())) {
00207         return reco::GenParticleRef(g.gens, it - g.genBarcodes->begin());
00208     } else {
00209         return reco::GenParticleRef();
00210     }
00211 }
00212 
00213 reco::GenParticle 
00214 PATGenCandsFromSimTracksProducer::makeGenParticle_(const SimTrack &tk, const edm::Ref<reco::GenParticleCollection> & mother, const GlobalContext &g) const {
00215     // Make up a GenParticleCandidate from the GEANT track info.
00216     int charge = static_cast<int>(tk.charge());
00217     Particle::LorentzVector p4 = tk.momentum();
00218     Particle::Point vtx; // = (0,0,0) by default
00219     if (!tk.noVertex()) vtx = g.simvtxs[tk.vertIndex()].position();
00220     GenParticle gp(charge, p4, vtx, tk.type(), setStatus_, true);
00221     if (mother.isNonnull()) gp.addMother(mother);
00222     return gp;
00223 }
00224 
00225 
00226 void PATGenCandsFromSimTracksProducer::produce(Event& event,
00227                                             const EventSetup& iSetup) {
00228 
00229   if (firstEvent_){
00230     if (!pdts_.empty()) {
00231       pdgIds_.clear();
00232       for (vector<PdtEntry>::iterator itp = pdts_.begin(), edp = pdts_.end(); itp != edp; ++itp) {
00233         itp->setup(iSetup); // decode string->pdgId and vice-versa                                                                                                
00234         pdgIds_.insert(std::abs(itp->pdgId()));
00235       }
00236       pdts_.clear();
00237     }
00238     if (!motherPdts_.empty()) {
00239       motherPdgIds_.clear();
00240       for (vector<PdtEntry>::iterator itp = motherPdts_.begin(), edp = motherPdts_.end(); itp != edp; ++itp) {
00241         itp->setup(iSetup); // decode string->pdgId and vice-versa                                                                                                
00242         motherPdgIds_.insert(std::abs(itp->pdgId()));
00243       }
00244       motherPdts_.clear();
00245     }
00246     firstEvent_ = false;
00247   }
00248 
00249   // Simulated tracks (i.e. GEANT particles).
00250   Handle<SimTrackContainer> simtracks;
00251   event.getByLabel(src_, simtracks);
00252 
00253   // Need to check that SimTrackContainer is sorted; otherwise, copy and sort :-(
00254   std::auto_ptr<SimTrackContainer> simtracksTmp;
00255   const SimTrackContainer * simtracksSorted = &* simtracks;
00256   if (makeMotherLink_ || writeAncestors_) {
00257       if (!__gnu_cxx::is_sorted(simtracks->begin(), simtracks->end(), LessById())) {
00258           simtracksTmp.reset(new SimTrackContainer(*simtracks));
00259           std::sort(simtracksTmp->begin(), simtracksTmp->end(), LessById());
00260           simtracksSorted = &* simtracksTmp;
00261       }
00262   }
00263 
00264   // Get the associated vertices
00265   Handle<SimVertexContainer> simvertices;
00266   event.getByLabel(src_, simvertices);
00267 
00268   // Get the GenParticles and barcodes, if needed to set mother links
00269   Handle<GenParticleCollection> gens;
00270   Handle<std::vector<int> > genBarcodes;
00271   bool barcodesAreSorted = true;
00272   if (makeMotherLink_) {
00273       event.getByLabel(genParticles_, gens);
00274       event.getByLabel(genParticles_, genBarcodes);
00275       if (gens->size() != genBarcodes->size()) throw cms::Exception("Corrupt data") << "Barcodes not of the same size as GenParticles!\n";
00276       barcodesAreSorted = __gnu_cxx::is_sorted(genBarcodes->begin(), genBarcodes->end());
00277   }    
00278 
00279 
00280   // make the output collection
00281   auto_ptr<GenParticleCollection> cands(new GenParticleCollection);
00282   edm::RefProd<GenParticleCollection> refprod = event.getRefBeforePut<GenParticleCollection>();
00283 
00284   GlobalContext globals(*simtracksSorted, *simvertices, gens, genBarcodes, barcodesAreSorted, *cands, refprod);
00285 
00286   for (SimTrackContainer::const_iterator isimtrk = simtracks->begin();
00287           isimtrk != simtracks->end(); ++isimtrk) {
00288 
00289       // Skip PYTHIA tracks.
00290       if (isimtrk->genpartIndex() != -1) continue; 
00291 
00292       // Maybe apply the PdgId filter
00293       if (!pdgIds_.empty()) { // if we have a filter on pdg ids
00294            if (pdgIds_.find(std::abs(isimtrk->type())) == pdgIds_.end()) continue;
00295       }
00296 
00297       GenParticle genp = makeGenParticle_(*isimtrk, Ref<GenParticleCollection>(), globals);
00298 
00299       // Maybe apply filter on the particle
00300       if (filter_.get() != 0) {
00301         if (!(*filter_)(genp)) continue;
00302       }
00303 
00304       if (!motherPdgIds_.empty()) {
00305            const SimTrack *motherSimTk = findGeantMother(*isimtrk, globals);
00306            if (motherSimTk == 0) continue;
00307            if (motherPdgIds_.find(std::abs(motherSimTk->type())) == motherPdgIds_.end()) continue;
00308       }
00309 
00310       if (makeMotherLink_ || writeAncestors_) {
00311           Ref<GenParticleCollection> motherRef;
00312           const SimTrack * mother = findGeantMother(*isimtrk, globals);
00313           if (mother != 0) motherRef = findRef(*mother, globals); 
00314           if (motherRef.isNonnull()) genp.addMother(motherRef);
00315       }
00316 
00317       cands->push_back(genp);
00318   }
00319 
00320   // Write to the Event, and get back a handle (which can be useful for debugging)
00321   edm::OrphanHandle<reco::GenParticleCollection> orphans = event.put(cands);
00322 
00323 #ifdef DEBUG_PATGenCandsFromSimTracksProducer
00324   std::cout << "Produced a list of " << orphans->size() << " genParticles." << std::endl;
00325   for (GenParticleCollection::const_iterator it = orphans->begin(), ed = orphans->end(); it != ed; ++it) {
00326       std::cout << "    ";
00327       std::cout << "GenParticle #" << (it - orphans->begin()) << ": pdgId " << it->pdgId()
00328                 << ", pt = " << it->pt() << ", eta = " << it->eta() << ", phi = " << it->phi() 
00329                 << ", rho = " << it->vertex().Rho() << ", z = " << it->vertex().Z() << std::endl;
00330       edm::Ref<GenParticleCollection> mom = it->motherRef();
00331       size_t depth = 2;
00332       while (mom.isNonnull()) {
00333           if (mom.id() == orphans.id()) {
00334               // I need to re-make the ref because they are not working until this module returns.
00335               mom = edm::Ref<GenParticleCollection>(orphans, mom.key());
00336           }
00337           for (size_t i = 0; i < depth; ++i) std::cout << "    ";
00338           std::cout << "GenParticleRef [" << mom.id() << "/" << mom.key() << "]: pdgId " << mom->pdgId() << ", status = " << mom->status() 
00339                     << ", pt = " << mom->pt() << ", eta = " << mom->eta() << ", phi = " << mom->phi() 
00340                     << ", rho = " << mom->vertex().Rho() << ", z = " << mom->vertex().Z() << std::endl;
00341           if (mom.id() != orphans.id()) break;
00342           if ((mom->motherRef().id() == mom.id()) && (mom->motherRef().key() == mom.key())) {
00343               throw cms::Exception("Corrupt Data") << "A particle is it's own mother.\n";
00344           }
00345           mom = mom->motherRef(); 
00346           depth++; 
00347       }
00348   }
00349   std::cout << std::endl;
00350 #endif
00351 
00352 }
00353 
00354 DEFINE_FWK_MODULE(PATGenCandsFromSimTracksProducer);