CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoJets/JetProducers/plugins/SubEventGenJetProducer.cc

Go to the documentation of this file.
00001 
00002 #include "FWCore/Framework/interface/MakerMacros.h"
00003 #include "RecoJets/JetProducers/plugins/SubEventGenJetProducer.h"
00004 #include "FWCore/Utilities/interface/Exception.h"
00005 #include "RecoJets/JetProducers/interface/JetSpecific.h"
00006 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00007 
00008 using namespace std;
00009 using namespace reco;
00010 using namespace edm;
00011 using namespace cms;
00012 
00013 namespace {
00014   const bool debug = false;
00015 
00016 }
00017 
00018 namespace {
00019    bool checkHydro(const reco::GenParticle * p){
00020       const Candidate* m1 = p->mother();
00021       while(m1){
00022          int pdg = abs(m1->pdgId());
00023          int st = m1->status();
00024          LogDebug("SubEventMothers")<<"Pdg ID : "<<pdg<<endl;
00025          if(st == 3 || pdg < 9 || pdg == 21){
00026             LogDebug("SubEventMothers")<<"Sub-Collision  Found! Pdg ID : "<<pdg<<endl;
00027             return false;
00028          }
00029          const Candidate* m = m1->mother();
00030          m1 = m;
00031          if(!m1) LogDebug("SubEventMothers")<<"No Mother, particle is : "<<pdg<<" with status "<<st<<endl;
00032       }
00033       //      return true;
00034       return true; // Debugging - to be changed
00035    }
00036 }
00037 
00038 SubEventGenJetProducer::SubEventGenJetProducer(edm::ParameterSet const& conf):
00039   VirtualJetProducer( conf )
00040 {
00041    //   mapSrc_ = conf.getParameter<edm::InputTag>( "srcMap");
00042    ignoreHydro_ = conf.getUntrackedParameter<bool>("ignoreHydro", true);
00043    produces<reco::BasicJetCollection>();
00044   // the subjet collections are set through the config file in the "jetCollInstanceName" field.
00045 }
00046 
00047 
00048 void SubEventGenJetProducer::inputTowers( ) 
00049 {
00050    std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
00051       inEnd = inputs_.end(), i = inBegin;
00052    for (; i != inEnd; ++i ) {
00053       reco::CandidatePtr input = inputs_[i - inBegin];
00054       if (isnan(input->pt()))           continue;
00055       if (input->et()    <inputEtMin_)  continue;
00056       if (input->energy()<inputEMin_)   continue;
00057       if (isAnomalousTower(input))      continue;
00058 
00059       edm::Ptr<reco::Candidate> p = inputs_[i - inBegin];
00060       const GenParticle * pref = dynamic_cast<const GenParticle *>(p.get());
00061       int subevent = pref->collisionId();
00062       LogDebug("SubEventContainers")<<"SubEvent is : "<<subevent<<endl;
00063       LogDebug("SubEventContainers")<<"SubSize is : "<<subInputs_.size()<<endl;
00064 
00065       if(subevent >= (int)subInputs_.size()){ 
00066          hydroTag_.resize(subevent+1, -1);
00067          subInputs_.resize(subevent+1);
00068       LogDebug("SubEventContainers")<<"SubSize is : "<<subInputs_.size()<<endl;
00069       LogDebug("SubEventContainers")<<"HydroTagSize is : "<<hydroTag_.size()<<endl;
00070       }
00071 
00072       LogDebug("SubEventContainers")<<"HydroTag is : "<<hydroTag_[subevent]<<endl;
00073       if(hydroTag_[subevent] != 0) hydroTag_[subevent] = (int)checkHydro(pref);
00074 
00075       subInputs_[subevent].push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
00076                                                 input->energy()));
00077 
00078       subInputs_[subevent].back().set_user_index(i - inBegin);
00079 
00080    }
00081 
00082 }
00083 
00084 void SubEventGenJetProducer::produce(edm::Event& iEvent,const edm::EventSetup& iSetup){
00085    LogDebug("VirtualJetProducer") << "Entered produce\n";
00086 
00087    fjJets_.clear();
00088    subInputs_.clear();
00089    nSubParticles_.clear();
00090    hydroTag_.clear();
00091    inputs_.clear();
00092 
00093    // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
00094    edm::Handle<reco::CandidateView> inputsHandle;
00095    iEvent.getByLabel(src_,inputsHandle);
00096    for (size_t i = 0; i < inputsHandle->size(); ++i) {
00097      inputs_.push_back(inputsHandle->ptrAt(i));
00098    }
00099    LogDebug("VirtualJetProducer") << "Got inputs\n";
00100 
00101    inputTowers();
00102    // Convert candidates to fastjet::PseudoJets.
00103    // Also correct to Primary Vertex. Will modify fjInputs_
00104    // and use inputs_
00105 
00107 
00108    std::auto_ptr<std::vector<GenJet> > jets(new std::vector<GenJet>() );
00109    subJets_ = jets.get();
00110 
00111    LogDebug("VirtualJetProducer") << "Inputted towers\n";
00112 
00113    size_t nsub = subInputs_.size();
00114 
00115    for(size_t isub = 0; isub < nsub; ++isub){
00116       if(ignoreHydro_ && hydroTag_[isub]) continue;
00117       fjJets_.clear();
00118       fjInputs_.clear();
00119       fjInputs_ = subInputs_[isub];
00120       runAlgorithm( iEvent, iSetup );
00121    }
00122 
00123    //Finalize
00124    LogDebug("SubEventJetProducer") << "Wrote jets\n";
00125 
00126    iEvent.put(jets);  
00127    return;
00128 }
00129 
00130 void SubEventGenJetProducer::runAlgorithm( edm::Event & iEvent, edm::EventSetup const& iSetup)
00131 {
00132    // run algorithm
00133    fjJets_.clear();
00134 
00135    fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
00136    fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
00137    
00138    using namespace reco;
00139 
00140    for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
00141       
00142       GenJet jet;
00143       const fastjet::PseudoJet& fjJet = fjJets_[ijet];
00144 
00145       std::vector<fastjet::PseudoJet> fjConstituents =
00146          sorted_by_pt(fjClusterSeq_->constituents(fjJet));
00147 
00148       std::vector<CandidatePtr> constituents =
00149          getConstituents(fjConstituents);
00150 
00151     double px=fjJet.px();
00152     double py=fjJet.py();
00153     double pz=fjJet.pz();
00154     double E=fjJet.E();
00155     double jetArea=0.0;
00156     double pu=0.;
00157 
00158     writeSpecific( jet,
00159                    Particle::LorentzVector(px, py, pz, E),
00160                    vertex_,
00161                    constituents, iSetup);
00162     
00163     jet.setJetArea (jetArea);
00164     jet.setPileup (pu);
00165     
00166     subJets_->push_back(jet);
00167    }   
00168 }
00169 
00170 DEFINE_FWK_MODULE(SubEventGenJetProducer);
00171 
00172