CMS 3D CMS Logo

Classes | Public Member Functions | Public Attributes | Private Attributes

PrimaryVertexProducer Class Reference

#include <RecoVertex/PrimaryVertexProducer/src/PrimaryVertexProducer.cc>

Inheritance diagram for PrimaryVertexProducer:
edm::EDProducer edm::ProducerBase edm::ProductRegistryHelper

List of all members.

Classes

struct  algo

Public Member Functions

edm::ParameterSet config () const
 PrimaryVertexProducer (const edm::ParameterSet &)
virtual void produce (edm::Event &, const edm::EventSetup &)
 ~PrimaryVertexProducer ()

Public Attributes

edm::InputTag beamSpotLabel
edm::InputTag trackLabel

Private Attributes

std::vector< algoalgorithms
bool fVerbose
edm::ParameterSet theConfig
TrackClusterizerInZtheTrackClusterizer
TrackFilterForPVFindingBasetheTrackFilter

Detailed Description

Description: steers tracker primary vertex reconstruction and storage

Implementation: <Notes on="" implementation>="">

Definition at line 56 of file PrimaryVertexProducer.h.


Constructor & Destructor Documentation

PrimaryVertexProducer::PrimaryVertexProducer ( const edm::ParameterSet conf) [explicit]

Definition at line 19 of file PrimaryVertexProducer.cc.

References algorithm(), algorithms, beamSpotLabel, edm::ParameterSet::exists(), PrimaryVertexProducer::algo::fitter, fVerbose, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), PrimaryVertexProducer::algo::label, PrimaryVertexProducer::algo::minNdof, theTrackClusterizer, theTrackFilter, trackLabel, PrimaryVertexProducer::algo::useBeamConstraint, and PrimaryVertexProducer::algo::vertexSelector.

  :theConfig(conf)
{

  fVerbose   = conf.getUntrackedParameter<bool>("verbose", false);
  trackLabel = conf.getParameter<edm::InputTag>("TrackLabel");
  beamSpotLabel = conf.getParameter<edm::InputTag>("beamSpotLabel");


  // select and configure the track selection
  std::string trackSelectionAlgorithm=conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
  if(trackSelectionAlgorithm=="filter"){
    theTrackFilter= new TrackFilterForPVFinding( conf.getParameter<edm::ParameterSet>("TkFilterParameters") );
  }else if (trackSelectionAlgorithm=="filterWithThreshold"){
    theTrackFilter= new HITrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
  }else{
    throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " + trackSelectionAlgorithm);  
  }


  // select and configure the track clusterizer
  std::string clusteringAlgorithm=conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
  if (clusteringAlgorithm=="gap"){
    theTrackClusterizer = new GapClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
  }else if(clusteringAlgorithm=="DA"){
    theTrackClusterizer = new DAClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
  }
  // provide the vectorized version of the clusterizer, if supported by the build
#ifdef __GXX_EXPERIMENTAL_CXX0X__
   else if(clusteringAlgorithm == "DA_vect") {
    theTrackClusterizer = new DAClusterizerInZ_vect(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
  }
#endif


  else{
    throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);  
  }


  // select and configure the vertex fitters
  if (conf.exists("vertexCollections")){
    std::vector<edm::ParameterSet> vertexCollections =conf.getParameter< std::vector<edm::ParameterSet> >("vertexCollections");

    for( std::vector< edm::ParameterSet >::const_iterator algoconf = vertexCollections.begin(); algoconf != vertexCollections.end(); algoconf++){
      
      algo algorithm;
      std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
      if (fitterAlgorithm=="KalmanVertexFitter") {
        algorithm.fitter= new KalmanVertexFitter();
      } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
        algorithm.fitter= new AdaptiveVertexFitter();
      } else {
        throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);  
      }
      algorithm.label = algoconf->getParameter<std::string>("label");
      algorithm.minNdof = algoconf->getParameter<double>("minNdof");
      algorithm.useBeamConstraint=algoconf->getParameter<bool>("useBeamConstraint");
      algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
      algorithms.push_back(algorithm);
      
      produces<reco::VertexCollection>(algorithm.label);
    }
  }else{
    edm::LogWarning("MisConfiguration")<<"this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";

    algo algorithm;
    std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
    if (fitterAlgorithm=="KalmanVertexFitter") {
      algorithm.fitter= new KalmanVertexFitter();
    } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
      algorithm.fitter= new AdaptiveVertexFitter();
    } else {
      throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);  
    }
    algorithm.label = "";
    algorithm.minNdof = conf.getParameter<double>("minNdof");
    algorithm.useBeamConstraint=conf.getParameter<bool>("useBeamConstraint");
    
    algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));

    algorithms.push_back(algorithm);
    produces<reco::VertexCollection>(algorithm.label);
  }
 

}
PrimaryVertexProducer::~PrimaryVertexProducer ( )

Definition at line 108 of file PrimaryVertexProducer.cc.

References algorithm(), algorithms, theTrackClusterizer, and theTrackFilter.

                                              {
  if (theTrackFilter) delete theTrackFilter;
  if (theTrackClusterizer) delete theTrackClusterizer;
  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
    if (algorithm->fitter) delete algorithm->fitter;
    if (algorithm->vertexSelector) delete algorithm->vertexSelector;
  }
}

Member Function Documentation

edm::ParameterSet PrimaryVertexProducer::config ( void  ) const [inline]

Definition at line 64 of file PrimaryVertexProducer.h.

References theConfig.

{ return theConfig; }
void PrimaryVertexProducer::produce ( edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDProducer.

Definition at line 119 of file PrimaryVertexProducer.cc.

References algorithm(), algorithms, SiPixelRawToDigiRegional_cfi::beamSpot, beamSpotLabel, TrackClusterizerInZ::clusterize(), gather_cfg::cout, GlobalErrorBase< T, ErrorWeightType >::cxx(), GlobalErrorBase< T, ErrorWeightType >::cyy(), GlobalErrorBase< T, ErrorWeightType >::czz(), TransientVertex::degreesOfFreedom(), VertexState::error(), fVerbose, edm::EventSetup::get(), edm::Event::getByLabel(), edm::HandleBase::isValid(), TransientVertex::isValid(), GlobalErrorBase< T, ErrorWeightType >::matrix(), TransientVertex::position(), reco::BeamSpot::position(), edm::Event::put(), query::result, reco::BeamSpot::rotatedCovariance3D(), TrackFilterForPVFindingBase::select(), python::multivaluedict::sort(), theTrackClusterizer, theTrackFilter, trackLabel, v, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

{

  // get the BeamSpot, it will alwys be needed, even when not used as a constraint
  reco::BeamSpot beamSpot;
  edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
  iEvent.getByLabel(beamSpotLabel,recoBeamSpotHandle);
  if (recoBeamSpotHandle.isValid()){
    beamSpot = *recoBeamSpotHandle;
  }else{
    edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
  }

  bool validBS = true;
  VertexState beamVertexState(beamSpot);
  if ( (beamVertexState.error().cxx() <= 0.) || 
       (beamVertexState.error().cyy() <= 0.) ||
       (beamVertexState.error().czz() <= 0.) ) {
    validBS = false;
    edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors "<<beamVertexState.error().matrix();
  }


  // get RECO tracks from the event
  // `tks` can be used as a ptr to a reco::TrackCollection
  edm::Handle<reco::TrackCollection> tks;
  iEvent.getByLabel(trackLabel, tks);


  // interface RECO tracks to vertex reconstruction
  edm::ESHandle<TransientTrackBuilder> theB;
  iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
  std::vector<reco::TransientTrack> t_tks = (*theB).build(tks, beamSpot);
  if(fVerbose) {std::cout << "RecoVertex/PrimaryVertexProducer"
                     << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
  }


  // select tracks
  std::vector<reco::TransientTrack> seltks = theTrackFilter->select( t_tks );


  // clusterize tracks in Z
  std::vector< std::vector<reco::TransientTrack> > clusters =  theTrackClusterizer->clusterize(seltks);
  if (fVerbose){std::cout <<  " clustering returned  "<< clusters.size() << " clusters  from " << seltks.size() << " selected tracks" <<std::endl;}


  // vertex fits
  for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){


    std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
    reco::VertexCollection vColl;


    std::vector<TransientVertex> pvs;
    for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
           = clusters.begin(); iclus != clusters.end(); iclus++) {


      TransientVertex v; 
      if( algorithm->useBeamConstraint && validBS &&((*iclus).size()>1) ){
        
        v = algorithm->fitter->vertex(*iclus, beamSpot);
        
      }else if( !(algorithm->useBeamConstraint) && ((*iclus).size()>1) ) {
      
        v = algorithm->fitter->vertex(*iclus); 
        
      }// else: no fit ==> v.isValid()=False


      if (fVerbose){
        if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " <<  v.position().z() << std::endl;
        else std::cout <<"Invalid fitted vertex\n";
      }

      if (v.isValid() 
            && (v.degreesOfFreedom()>=algorithm->minNdof) 
          && (!validBS || (*(algorithm->vertexSelector))(v,beamVertexState))
          ) pvs.push_back(v);
    }// end of cluster loop

    if(fVerbose){
      std::cout << "PrimaryVertexProducerAlgorithm::vertices  candidates =" << pvs.size() << std::endl;
    }



    

    // sort vertices by pt**2  vertex (aka signal vertex tagging)
    if(pvs.size()>1){
      sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
    }



    // convert transient vertices returned by the theAlgo to (reco) vertices
    for (std::vector<TransientVertex>::const_iterator iv = pvs.begin();
         iv != pvs.end(); iv++) {
      reco::Vertex v = *iv;
      vColl.push_back(v);
    }

    if (vColl.empty()) {
      GlobalError bse(beamSpot.rotatedCovariance3D());
      if ( (bse.cxx() <= 0.) || 
           (bse.cyy() <= 0.) ||
           (bse.czz() <= 0.) ) {
        AlgebraicSymMatrix33 we;
        we(0,0)=10000; we(1,1)=10000; we(2,2)=10000;
        vColl.push_back(reco::Vertex(beamSpot.position(), we,0.,0.,0));
        if(fVerbose){
          std::cout <<"RecoVertex/PrimaryVertexProducer: "
                    << "Beamspot with invalid errors "<<bse.matrix()<<std::endl;
          std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
        }
      } else {
        vColl.push_back(reco::Vertex(beamSpot.position(), 
                                     beamSpot.rotatedCovariance3D(),0.,0.,0));
        if(fVerbose){
          std::cout <<"RecoVertex/PrimaryVertexProducer: "
                    << " will put Vertex derived from BeamSpot into Event.\n";
        }
      }
    }

    if(fVerbose){
      int ivtx=0;
      for(reco::VertexCollection::const_iterator v=vColl.begin(); 
          v!=vColl.end(); ++v){
        std::cout << "recvtx "<< ivtx++ 
                  << "#trk " << std::setw(3) << v->tracksSize()
                  << " chi2 " << std::setw(4) << v->chi2() 
                  << " ndof " << std::setw(3) << v->ndof() 
                  << " x "  << std::setw(6) << v->position().x() 
                  << " dx " << std::setw(6) << v->xError()
                  << " y "  << std::setw(6) << v->position().y() 
                  << " dy " << std::setw(6) << v->yError()
                  << " z "  << std::setw(6) << v->position().z() 
                  << " dz " << std::setw(6) << v->zError()
                  << std::endl;
      }
    }

  
    *result = vColl;
    iEvent.put(result, algorithm->label); 
  }
  
}

Member Data Documentation

std::vector< algo > PrimaryVertexProducer::algorithms [private]

Definition at line 82 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), produce(), and ~PrimaryVertexProducer().

Definition at line 66 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

Definition at line 85 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().

Definition at line 84 of file PrimaryVertexProducer.h.

Referenced by config().

Definition at line 71 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), produce(), and ~PrimaryVertexProducer().

Definition at line 70 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), produce(), and ~PrimaryVertexProducer().

Definition at line 65 of file PrimaryVertexProducer.h.

Referenced by PrimaryVertexProducer(), and produce().