CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Private Attributes

FastjetJetProducer Class Reference

#include <FastjetJetProducer.h>

Inheritance diagram for FastjetJetProducer:
VirtualJetProducer edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Public Member Functions

 FastjetJetProducer (const edm::ParameterSet &iConfig)
virtual void produce (edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual ~FastjetJetProducer ()

Protected Member Functions

virtual void produceTrackJets (edm::Event &iEvent, const edm::EventSetup &iSetup)
virtual void runAlgorithm (edm::Event &iEvent, const edm::EventSetup &iSetup)

Private Attributes

float dxyTrVtxMax_
float dzTrVtxMax_
float maxVtxZ_
int minVtxNdof_
bool useOnlyOnePV_
bool useOnlyVertexTracks_

Detailed Description

Definition at line 8 of file FastjetJetProducer.h.


Constructor & Destructor Documentation

FastjetJetProducer::FastjetJetProducer ( const edm::ParameterSet iConfig) [explicit]

Definition at line 56 of file FastjetJetProducer.cc.

References dxyTrVtxMax_, dzTrVtxMax_, edm::ParameterSet::exists(), edm::ParameterSet::getParameter(), maxVtxZ_, minVtxNdof_, useOnlyOnePV_, and useOnlyVertexTracks_.

  : VirtualJetProducer( iConfig )
{
  
  if ( iConfig.exists("UseOnlyVertexTracks") )
    useOnlyVertexTracks_ = iConfig.getParameter<bool>("UseOnlyVertexTracks");
  else 
    useOnlyVertexTracks_ = false;
  
  if ( iConfig.exists("UseOnlyOnePV") )
    useOnlyOnePV_        = iConfig.getParameter<bool>("UseOnlyOnePV");
  else
    useOnlyOnePV_ = false;

  if ( iConfig.exists("DzTrVtxMax") )
    dzTrVtxMax_          = iConfig.getParameter<double>("DzTrVtxMax");
  else
    dzTrVtxMax_ = 999999.;
  if ( iConfig.exists("DxyTrVtxMax") )
    dxyTrVtxMax_          = iConfig.getParameter<double>("DxyTrVtxMax");
  else
    dxyTrVtxMax_ = 999999.;
  if ( iConfig.exists("MinVtxNdof") )
    minVtxNdof_ = iConfig.getParameter<int>("MinVtxNdof");
  else
    minVtxNdof_ = 5;
  if ( iConfig.exists("MaxVtxZ") )
    maxVtxZ_ = iConfig.getParameter<double>("MaxVtxZ");
  else
    maxVtxZ_ = 15;

}
FastjetJetProducer::~FastjetJetProducer ( ) [virtual]

Definition at line 91 of file FastjetJetProducer.cc.

{
} 

Member Function Documentation

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

Reimplemented from VirtualJetProducer.

Definition at line 100 of file FastjetJetProducer.cc.

References VirtualJetProducer::jetTypeE, VirtualJetProducer::makeTrackJet(), and produceTrackJets().

{

  // for everything but track jets
  if (!makeTrackJet(jetTypeE)) {
 
     // use the default production from one collection
     VirtualJetProducer::produce( iEvent, iSetup );

  } else { // produce trackjets from tracks grouped per primary vertex

    produceTrackJets(iEvent, iSetup);
  
  }

}
void FastjetJetProducer::produceTrackJets ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [protected, virtual]

Definition at line 118 of file FastjetJetProducer.cc.

References dxyTrVtxMax_, dzTrVtxMax_, VirtualJetProducer::fjClusterSeq_, VirtualJetProducer::fjInputs_, VirtualJetProducer::fjJets_, newFWLiteAna::found, edm::Event::getByLabel(), VirtualJetProducer::getConstituents(), i, VirtualJetProducer::inputs_, VirtualJetProducer::inputTowers(), metsig::jet, analyzePatCleaning_cfg::jets, LogDebug, maxVtxZ_, minVtxNdof_, edm::Event::put(), runAlgorithm(), reco::Jet::setJetArea(), reco::Jet::setPileup(), reco::TrackJet::setPrimaryVertex(), reco::LeafCandidate::setVertex(), VirtualJetProducer::src_, VirtualJetProducer::srcPVs_, useOnlyOnePV_, useOnlyVertexTracks_, VirtualJetProducer::vertex_, and reco::writeSpecific().

Referenced by produce().

{

    // read in the track candidates
    edm::Handle<edm::View<reco::RecoChargedRefCandidate> > inputsHandle;
    iEvent.getByLabel(src_, inputsHandle);
    // make collection with pointers so we can play around with it
    std::vector<edm::Ptr<reco::RecoChargedRefCandidate> > allInputs;
    std::vector<edm::Ptr<reco::Candidate> > origInputs;
    for (size_t i = 0; i < inputsHandle->size(); ++i) {
      allInputs.push_back(inputsHandle->ptrAt(i));
      origInputs.push_back(inputsHandle->ptrAt(i));
    }

    // read in the PV collection
    edm::Handle<reco::VertexCollection> pvCollection;
    iEvent.getByLabel(srcPVs_, pvCollection);
    // define the overall output jet container
    std::auto_ptr<std::vector<reco::TrackJet> > jets(new std::vector<reco::TrackJet>() );

    // loop over the good vertices, clustering for each vertex separately
    for (reco::VertexCollection::const_iterator itVtx = pvCollection->begin(); itVtx != pvCollection->end(); ++itVtx) {
      if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;

      // clear the intermediate containers
      inputs_.clear();
      fjInputs_.clear();
      fjJets_.clear();

      // if only vertex-associated tracks should be used
      if (useOnlyVertexTracks_) {
        // loop over the tracks associated to the vertex
        for (reco::Vertex::trackRef_iterator itTr = itVtx->tracks_begin(); itTr != itVtx->tracks_end(); ++itTr) {
          // whether a match was found in the track candidate input
          bool found = false;
          // loop over input track candidates
          for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
            // match the input track candidate to the track from the vertex
            reco::TrackRef trref(itTr->castTo<reco::TrackRef>());
            // check if the tracks match
            if ((*itIn)->track() == trref) {
              found = true;
              // add this track candidate to the input for clustering
              inputs_.push_back(*itIn);
              // erase the track candidate from the total list of input, so we don't reuse it later
              allInputs.erase(itIn);
              // found the candidate track corresponding to the vertex track, so stop the loop
              break;
            } // end if match found
          } // end loop over input tracks
          // give an info message in case no match is found (can happen if candidates are subset of tracks used for clustering)
          if (!found) edm::LogInfo("FastjetTrackJetProducer") << "Ignoring a track at vertex which is not in input track collection!";
        } // end loop over tracks associated to vertex
      // if all inpt track candidates should be used
      } else {
        // loop over input track candidates
        for (std::vector<edm::Ptr<reco::RecoChargedRefCandidate> >::iterator itIn = allInputs.begin(); itIn != allInputs.end(); ++itIn) {
          // check if the track is close enough to the vertex
          float dz = (*itIn)->track()->dz(itVtx->position());
          float dxy = (*itIn)->track()->dxy(itVtx->position());
          if (fabs(dz) > dzTrVtxMax_) continue;
          if (fabs(dxy) > dxyTrVtxMax_) continue;
          bool closervtx = false;
          // now loop over the good vertices a second time
          for (reco::VertexCollection::const_iterator itVtx2 = pvCollection->begin(); itVtx2 != pvCollection->end(); ++itVtx2) {
            if (itVtx->isFake() || itVtx->ndof() < minVtxNdof_ || fabs(itVtx->z()) > maxVtxZ_) continue;
            // and check this track is closer to any other vertex (if more than 1 vertex considered)
            if (!useOnlyOnePV_ &&
                itVtx != itVtx2 &&
                fabs((*itIn)->track()->dz(itVtx2->position())) < fabs(dz)) {
              closervtx = true;
              break; // 1 closer vertex makes the track already not matched, so break
            }
          }
          // don't add this track if another vertex is found closer
          if (closervtx) continue;
          // add this track candidate to the input for clustering
          inputs_.push_back(*itIn);
          // erase the track candidate from the total list of input, so we don't reuse it later
          allInputs.erase(itIn);
          // take a step back in the loop since we just erased
          --itIn;
        }
      }

      // convert candidates in inputs_ to fastjet::PseudoJets in fjInputs_
      fjInputs_.reserve(inputs_.size());
      inputTowers();
      LogDebug("FastjetTrackJetProducer") << "Inputted towers\n";

      // run algorithm, using fjInputs_, modifying fjJets_ and allocating fjClusterSeq_
      runAlgorithm(iEvent, iSetup);
      LogDebug("FastjetTrackJetProducer") << "Ran algorithm\n";

      // convert our jets and add to the overall jet vector
      for (unsigned int ijet=0;ijet<fjJets_.size();++ijet) {
        // get the constituents from fastjet
        std::vector<fastjet::PseudoJet> fjConstituents = sorted_by_pt(fjClusterSeq_->constituents(fjJets_[ijet]));
        // convert them to CandidatePtr vector
        std::vector<reco::CandidatePtr> constituents = getConstituents(fjConstituents);
        // fill the trackjet
        reco::TrackJet jet;
        // write the specifics to the jet (simultaneously sets 4-vector, vertex).
        writeSpecific( jet,
                       reco::Particle::LorentzVector(fjJets_[ijet].px(), fjJets_[ijet].py(), fjJets_[ijet].pz(), fjJets_[ijet].E()),
                       vertex_, constituents, iSetup);
        jet.setJetArea(0);
        jet.setPileup(0);
        jet.setPrimaryVertex(edm::Ref<reco::VertexCollection>(pvCollection, (int) (itVtx-pvCollection->begin())));
        jet.setVertex(itVtx->position());
        jets->push_back(jet);
      }

      if (useOnlyOnePV_) break; // stop vertex loop if only one vertex asked for
    } // end loop over vertices

    // put the jets in the collection
    LogDebug("FastjetTrackJetProducer") << "Put " << jets->size() << " jets in the event.\n";
    iEvent.put(jets);

}
void FastjetJetProducer::runAlgorithm ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [protected, virtual]

Implements VirtualJetProducer.

Definition at line 242 of file FastjetJetProducer.cc.

References VirtualJetProducer::doAreaFastjet_, VirtualJetProducer::doRhoFastjet_, VirtualJetProducer::fjActiveArea_, VirtualJetProducer::fjClusterSeq_, VirtualJetProducer::fjInputs_, VirtualJetProducer::fjJetDefinition_, VirtualJetProducer::fjJets_, VirtualJetProducer::jetPtMin_, and VirtualJetProducer::voronoiRfact_.

Referenced by produceTrackJets().

{
  // run algorithm
  /*
  fjInputs_.clear();
  double px, py , pz, E;
  string line;
  std::ifstream fin("dump3.txt");
  while (getline(fin, line)){
    if (line == "#END") break;
    if (line.substr(0,1) == "#") {continue;}
    istringstream istr(line);
    istr >> px >> py >> pz >> E;
    // create a fastjet::PseudoJet with these components and put it onto
    // back of the input_particles vector
    fastjet::PseudoJet j(px,py,pz,E);
    //if ( fabs(j.rap()) < inputEtaMax )
    fjInputs_.push_back(fastjet::PseudoJet(px,py,pz,E)); 
  }
  fin.close();
  */
  if ( !doAreaFastjet_ && !doRhoFastjet_) {
    fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequence( fjInputs_, *fjJetDefinition_ ) );
  } else if (voronoiRfact_ <= 0) {
    fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceArea( fjInputs_, *fjJetDefinition_ , *fjActiveArea_ ) );
  } else {
    fjClusterSeq_ = ClusterSequencePtr( new fastjet::ClusterSequenceVoronoiArea( fjInputs_, *fjJetDefinition_ , fastjet::VoronoiAreaSpec(voronoiRfact_) ) );
  }
  /*
    std::cout << "using " << fjInputs_.size() << " particles.\n";
    std::cout.precision(16);//std::numeric_limits< double >::digits10);
    std::cout << "dumping input px py pz E:\n";
    for(std::vector<fastjet::PseudoJet>::const_iterator i = fjInputs_.begin() ; i != fjInputs_.end(); ++i) {
    std::cout << i->px() << " " << i->py() << " " << i->pz() << " " << i->E() << '\n';
    }
  */
  fjJets_ = fastjet::sorted_by_pt(fjClusterSeq_->inclusive_jets(jetPtMin_));
}

Member Data Documentation

Definition at line 36 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

Definition at line 35 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

Definition at line 38 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

Definition at line 37 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

Definition at line 34 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().

Definition at line 33 of file FastjetJetProducer.h.

Referenced by FastjetJetProducer(), and produceTrackJets().