CMS 3D CMS Logo

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