CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoVertex/PrimaryVertexProducer/src/PrimaryVertexProducerAlgorithm.cc

Go to the documentation of this file.
00001 #include "RecoVertex/PrimaryVertexProducer/interface/PrimaryVertexProducerAlgorithm.h"
00002 #include "RecoVertex/PrimaryVertexProducer/interface/VertexHigherPtSquared.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00006 #include "RecoVertex/PrimaryVertexProducer/interface/TrackFilterForPVFinding.h"
00007 #include "RecoVertex/PrimaryVertexProducer/interface/HITrackFilterForPVFinding.h"
00008 #include "RecoVertex/PrimaryVertexProducer/interface/GapClusterizerInZ.h"
00009 #include "RecoVertex/PrimaryVertexProducer/interface/DAClusterizerInZ.h"
00010 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
00011 #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
00012 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
00013 #include "RecoVertex/VertexPrimitives/interface/VertexException.h"
00014 #include <algorithm>
00015 
00016 using namespace reco;
00017 
00018 //
00019 // constructors and destructor
00020 //
00021 PrimaryVertexProducerAlgorithm::PrimaryVertexProducerAlgorithm(const edm::ParameterSet& conf)
00022   // extract relevant parts of config for components
00023   : theConfig(conf), 
00024     theVertexSelector(VertexDistanceXY(), 
00025                       conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"))
00026 {
00027   edm::LogInfo("PVDebugInfo") 
00028     << "PVSelParameters::maxDistanceToBeam = " 
00029     << conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam") << "\n";
00030 
00031 
00032   fUseBeamConstraint = conf.getParameter<bool>("useBeamConstraint");
00033   fVerbose           = conf.getUntrackedParameter<bool>("verbose", false);
00034   fMinNdof           = conf.getParameter<double>("minNdof");
00035   fFailsafe          = true; //conf.getUntrackedParameter<bool>("failsafe",true);
00036 
00037 
00038   // select and configure the track selection
00039   std::string trackSelectionAlgorithm=conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
00040   if(trackSelectionAlgorithm=="filter"){
00041     theTrackFilter= new TrackFilterForPVFinding( conf.getParameter<edm::ParameterSet>("TkFilterParameters") );
00042   }else if (trackSelectionAlgorithm=="filterWithThreshold"){
00043     theTrackFilter= new HITrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
00044   }else{
00045     throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " + trackSelectionAlgorithm);  
00046   }
00047 
00048 
00049   // select and configure the track clusterizer
00050   std::string clusteringAlgorithm=conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
00051   if (clusteringAlgorithm=="gap"){
00052     theTrackClusterizer = new GapClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
00053   }else if(clusteringAlgorithm=="DA"){
00054     theTrackClusterizer = new DAClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
00055   }else{
00056     throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);  
00057   }
00058 
00059   // select and configure the vertex fitter
00060   std::string algorithm = conf.getParameter<std::string>("algorithm");
00061   fapply_finder = false;
00062   if (algorithm == "TrimmedKalmanFinder") {
00063     fapply_finder = true;
00064     theFinder.setParameters(conf.getParameter<edm::ParameterSet>("VtxFinderParameters"));
00065   } else if (algorithm=="KalmanVertexFitter") {
00066     theFitter=new KalmanVertexFitter();
00067   } else if( algorithm=="AdaptiveVertexFitter") {
00068     theFitter=new AdaptiveVertexFitter();
00069   } else {
00070     throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + algorithm);  
00071   }
00072 
00073   edm::LogInfo("PVDebugInfo") 
00074     << "Using " << algorithm << "\n";
00075   edm::LogInfo("PVDebugInfo") 
00076     << "beam-constraint  " << fUseBeamConstraint << "\n"; 
00077 
00078   edm::LogInfo("PVDebugInfo") 
00079     << "PV producer algorithm initialization: done" << "\n";
00080 
00081 }
00082 
00083 
00084 PrimaryVertexProducerAlgorithm::~PrimaryVertexProducerAlgorithm() 
00085 {
00086   if (theFitter) delete theFitter;
00087   if (theTrackFilter) delete theTrackFilter;
00088   if (theTrackClusterizer) delete theTrackClusterizer;
00089 }
00090 
00091 
00092 //
00093 // member functions
00094 //
00095 
00096 // obsolete method, unfortunately required throgh inheritance from  VertexReconstructor
00097 std::vector<TransientVertex> 
00098 PrimaryVertexProducerAlgorithm::vertices(const std::vector<reco::TransientTrack> & tracks) const
00099 {
00100 
00101    throw VertexException("PrimaryVertexProducerAlgorithm: cannot make a Primary Vertex without a beam spot constraint " );
00102 
00103   /*  std::cout<< "PrimaryVertexProducer::vertices> Obsolete function, using dummy beamspot " << std::endl;
00104     reco::BeamSpot dummyBeamSpot;
00105     dummyBeamSpot.dummy();
00106     return vertices(tracks,dummyBeamSpot); */
00107    return std::vector<TransientVertex>();
00108 }
00109 
00110 
00111 std::vector<TransientVertex> 
00112 PrimaryVertexProducerAlgorithm::vertices(const std::vector<reco::TransientTrack> & tracks,
00113                                          const reco::BeamSpot & beamSpot) const
00114 {
00115   bool validBS = true;
00116   VertexState beamVertexState(beamSpot);
00117   if ( (beamVertexState.error().cxx() <= 0.) || 
00118        (beamVertexState.error().cyy() <= 0.) ||
00119        (beamVertexState.error().czz() <= 0.) ) {
00120     validBS = false;
00121     edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors "<<beamVertexState.error().matrix();
00122   }
00123 
00124   if ( fapply_finder) {
00125         return theFinder.vertices( tracks );
00126   }
00127   std::vector<TransientVertex> pvs;
00128 
00129 
00130   // select tracks
00131   std::vector<TransientTrack> seltks = theTrackFilter->select( tracks );
00132 
00133 
00134   // clusterize tracks in Z
00135   std::vector< std::vector<reco::TransientTrack> > clusters =  theTrackClusterizer->clusterize(seltks);
00136   if (fVerbose){std::cout <<  " clustering returned  "<< clusters.size() << " clusters  from " << seltks.size() << " selected tracks" <<std::endl;}
00137 
00138 
00139   // look for primary vertices in each cluster
00140   std::vector<TransientVertex> pvCand;
00141   int nclu=0;
00142   for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
00143          = clusters.begin(); iclus != clusters.end(); iclus++) {
00144 
00145 
00146     TransientVertex v;
00147     if( fUseBeamConstraint && validBS &&((*iclus).size()>1) ){
00148       if (fVerbose){std::cout <<  " constrained fit with "<< (*iclus).size() << " tracks"  <<std::endl;}
00149       v = theFitter->vertex(*iclus, beamSpot);
00150       if (v.isValid() && (v.degreesOfFreedom()>=fMinNdof)) pvCand.push_back(v);
00151 
00152       if (fVerbose){
00153         if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " <<  v.position().z() << std::endl;
00154         else std::cout <<"Invalid fitted vertex\n";
00155       }
00156 
00157     }else if((*iclus).size()>1){
00158       if (fVerbose){std::cout <<  " unconstrained fit with "<< (*iclus).size() << " tracks"  << std::endl;}
00159 
00160       v = theFitter->vertex(*iclus); 
00161       if (v.isValid() && (v.degreesOfFreedom()>=fMinNdof)) pvCand.push_back(v);
00162 
00163       if (fVerbose){
00164         if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " <<  v.position().z() << std::endl;
00165         else std::cout <<"Invalid fitted vertex\n";
00166       }
00167 
00168     }
00169 
00170     nclu++;
00171 
00172   }// end of cluster loop
00173 
00174   if(fVerbose){
00175     std::cout << "PrimaryVertexProducerAlgorithm::vertices  candidates =" << pvCand.size() << std::endl;
00176   }
00177 
00178 
00179 
00180   // select vertices compatible with beam
00181   int npv=0;
00182   for (std::vector<TransientVertex>::const_iterator ipv = pvCand.begin();
00183        ipv != pvCand.end(); ipv++) {
00184     if(fVerbose){
00185       std::cout << "PrimaryVertexProducerAlgorithm::vertices cand " << npv++ << " sel=" <<
00186         (validBS && theVertexSelector(*ipv,beamVertexState)) << "   z="  << ipv->position().z() << std::endl;
00187     }
00188     if (!validBS || theVertexSelector(*ipv,beamVertexState)) pvs.push_back(*ipv);
00189   }
00190 
00191 
00192   if(pvs.size()>0){
00193 
00194     // sort vertices by pt**2  vertex (aka signal vertex tagging)
00195     sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
00196 
00197   }else{
00198 
00199     if (    fFailsafe 
00200             && (seltks.size()>1) 
00201             && ( (clusters.size()!=1)  ||  ( (clusters.size()==1) && (clusters.begin()->size()<seltks.size())) )
00202         )
00203       { 
00204         // if no vertex was found, try fitting all selected tracks, unless this has already been tried
00205         // in low/no pile-up situations with low multiplicity vertices, this can recover vertices lost in clustering
00206         // with very small zSep. Only makes sense when used with a robust fitter, like the AdaptiveVertexFitter
00207         TransientVertex v;
00208         if( fUseBeamConstraint && validBS ){
00209           v = theFitter->vertex(seltks, beamSpot);
00210         }else{
00211           v = theFitter->vertex(seltks); 
00212         }
00213         if (fVerbose){ std::cout << "PrimaryVertexProducerAlgorithm: failsafe vertex "
00214                             <<"  tracks="  << seltks.size()
00215                             <<"  valid()=" << v.isValid() << " ndof()=" << v.degreesOfFreedom() 
00216                             <<"  selected="<< theVertexSelector(v,beamVertexState) << std::endl; }
00217         if ( v.isValid() && (v.degreesOfFreedom()>=fMinNdof) && (theVertexSelector(v,beamVertexState)) )  pvs.push_back(v);
00218       }
00219     
00220   }
00221 
00222 
00223   return pvs;
00224   
00225 }