CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Protected Attributes

cms::SubEventGenJetProducer Class Reference

#include <SubEventGenJetProducer.h>

Inheritance diagram for cms::SubEventGenJetProducer:
VirtualJetProducer edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

void produce (edm::Event &, const edm::EventSetup &)
void runAlgorithm (edm::Event &, const edm::EventSetup &)
 SubEventGenJetProducer (const edm::ParameterSet &ps)
virtual ~SubEventGenJetProducer ()

Protected Member Functions

virtual void inputTowers ()

Protected Attributes

std::vector< int > hydroTag_
bool ignoreHydro_
std::vector< int > nSubParticles_
std::vector< std::vector
< fastjet::PseudoJet > > 
subInputs_
std::vector< reco::GenJet > * subJets_

Detailed Description

Definition at line 19 of file SubEventGenJetProducer.h.


Constructor & Destructor Documentation

SubEventGenJetProducer::SubEventGenJetProducer ( const edm::ParameterSet ps)

Definition at line 38 of file SubEventGenJetProducer.cc.

References edm::ParameterSet::getUntrackedParameter(), and ignoreHydro_.

                                                                         :
  VirtualJetProducer( conf )
{
   //   mapSrc_ = conf.getParameter<edm::InputTag>( "srcMap");
   ignoreHydro_ = conf.getUntrackedParameter<bool>("ignoreHydro", true);
   produces<reco::BasicJetCollection>();
  // the subjet collections are set through the config file in the "jetCollInstanceName" field.
}
virtual cms::SubEventGenJetProducer::~SubEventGenJetProducer ( ) [inline, virtual]

Definition at line 24 of file SubEventGenJetProducer.h.

{}

Member Function Documentation

void SubEventGenJetProducer::inputTowers ( ) [protected, virtual]

Reimplemented from VirtualJetProducer.

Definition at line 48 of file SubEventGenJetProducer.cc.

References reco::GenParticle::collisionId(), edm::Ptr< T >::get(), hydroTag_, i, LaserDQM_cfg::input, VirtualJetProducer::inputEMin_, VirtualJetProducer::inputEtMin_, VirtualJetProducer::inputs_, VirtualJetProducer::isAnomalousTower(), edm::detail::isnan(), LogDebug, and subInputs_.

Referenced by produce().

{
   std::vector<edm::Ptr<reco::Candidate> >::const_iterator inBegin = inputs_.begin(),
      inEnd = inputs_.end(), i = inBegin;
   for (; i != inEnd; ++i ) {
      reco::CandidatePtr input = inputs_[i - inBegin];
      if (std::isnan(input->pt()))           continue;
      if (input->et()    <inputEtMin_)  continue;
      if (input->energy()<inputEMin_)   continue;
      if (isAnomalousTower(input))      continue;

      edm::Ptr<reco::Candidate> p = inputs_[i - inBegin];
      const GenParticle * pref = dynamic_cast<const GenParticle *>(p.get());
      int subevent = pref->collisionId();
      LogDebug("SubEventContainers")<<"SubEvent is : "<<subevent<<endl;
      LogDebug("SubEventContainers")<<"SubSize is : "<<subInputs_.size()<<endl;

      if(subevent >= (int)subInputs_.size()){ 
         hydroTag_.resize(subevent+1, -1);
         subInputs_.resize(subevent+1);
      LogDebug("SubEventContainers")<<"SubSize is : "<<subInputs_.size()<<endl;
      LogDebug("SubEventContainers")<<"HydroTagSize is : "<<hydroTag_.size()<<endl;
      }

      LogDebug("SubEventContainers")<<"HydroTag is : "<<hydroTag_[subevent]<<endl;
      if(hydroTag_[subevent] != 0) hydroTag_[subevent] = (int)checkHydro(pref);

      subInputs_[subevent].push_back(fastjet::PseudoJet(input->px(),input->py(),input->pz(),
                                                input->energy()));

      subInputs_[subevent].back().set_user_index(i - inBegin);

   }

}
void SubEventGenJetProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Reimplemented from VirtualJetProducer.

Definition at line 84 of file SubEventGenJetProducer.cc.

References VirtualJetProducer::fjInputs_, VirtualJetProducer::fjJets_, edm::Event::getByLabel(), hydroTag_, i, ignoreHydro_, VirtualJetProducer::inputs_, inputTowers(), fwrapper::jets, LogDebug, nSubParticles_, edm::Event::put(), runAlgorithm(), VirtualJetProducer::src_, subInputs_, and subJets_.

                                                                                {
   LogDebug("VirtualJetProducer") << "Entered produce\n";

   fjJets_.clear();
   subInputs_.clear();
   nSubParticles_.clear();
   hydroTag_.clear();
   inputs_.clear();

   // get inputs and convert them to the fastjet format (fastjet::PeudoJet)
   edm::Handle<reco::CandidateView> inputsHandle;
   iEvent.getByLabel(src_,inputsHandle);
   for (size_t i = 0; i < inputsHandle->size(); ++i) {
     inputs_.push_back(inputsHandle->ptrAt(i));
   }
   LogDebug("VirtualJetProducer") << "Got inputs\n";

   inputTowers();
   // Convert candidates to fastjet::PseudoJets.
   // Also correct to Primary Vertex. Will modify fjInputs_
   // and use inputs_


   std::auto_ptr<std::vector<GenJet> > jets(new std::vector<GenJet>() );
   subJets_ = jets.get();

   LogDebug("VirtualJetProducer") << "Inputted towers\n";

   size_t nsub = subInputs_.size();

   for(size_t isub = 0; isub < nsub; ++isub){
      if(ignoreHydro_ && hydroTag_[isub]) continue;
      fjJets_.clear();
      fjInputs_.clear();
      fjInputs_ = subInputs_[isub];
      runAlgorithm( iEvent, iSetup );
   }

   //Finalize
   LogDebug("SubEventJetProducer") << "Wrote jets\n";

   iEvent.put(jets);  
   return;
}
void SubEventGenJetProducer::runAlgorithm ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements VirtualJetProducer.

Definition at line 130 of file SubEventGenJetProducer.cc.

References VirtualJetProducer::fjClusterSeq_, VirtualJetProducer::fjInputs_, VirtualJetProducer::fjJetDefinition_, VirtualJetProducer::fjJets_, VirtualJetProducer::getConstituents(), metsig::jet, VirtualJetProducer::jetPtMin_, dt_dqm_sourceclient_common_cff::reco, reco::Jet::setJetArea(), reco::Jet::setPileup(), subJets_, VirtualJetProducer::vertex_, and reco::writeSpecific().

Referenced by produce().

{
   // run algorithm
   fjJets_.clear();

   fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
   fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
   
   using namespace reco;

   for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
      
      GenJet jet;
      const fastjet::PseudoJet& fjJet = fjJets_[ijet];

      std::vector<fastjet::PseudoJet> fjConstituents =
         sorted_by_pt(fjClusterSeq_->constituents(fjJet));

      std::vector<CandidatePtr> constituents =
         getConstituents(fjConstituents);

    double px=fjJet.px();
    double py=fjJet.py();
    double pz=fjJet.pz();
    double E=fjJet.E();
    double jetArea=0.0;
    double pu=0.;

    writeSpecific( jet,
                   Particle::LorentzVector(px, py, pz, E),
                   vertex_,
                   constituents, iSetup);
    
    jet.setJetArea (jetArea);
    jet.setPileup (pu);
    
    subJets_->push_back(jet);
   }   
}

Member Data Documentation

std::vector<int> cms::SubEventGenJetProducer::hydroTag_ [protected]

Definition at line 31 of file SubEventGenJetProducer.h.

Referenced by inputTowers(), and produce().

Definition at line 33 of file SubEventGenJetProducer.h.

Referenced by produce(), and SubEventGenJetProducer().

std::vector<int> cms::SubEventGenJetProducer::nSubParticles_ [protected]

Definition at line 32 of file SubEventGenJetProducer.h.

Referenced by produce().

std::vector<std::vector<fastjet::PseudoJet> > cms::SubEventGenJetProducer::subInputs_ [protected]

Definition at line 29 of file SubEventGenJetProducer.h.

Referenced by inputTowers(), and produce().

Definition at line 30 of file SubEventGenJetProducer.h.

Referenced by produce(), and runAlgorithm().