CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoVertex/PrimaryVertexProducer/plugins/PrimaryVertexProducer.cc

Go to the documentation of this file.
00001 #include "RecoVertex/PrimaryVertexProducer/interface/PrimaryVertexProducer.h"
00002 
00003 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00004 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00005 #include "DataFormats/Common/interface/Handle.h"
00006 #include "FWCore/Framework/interface/MakerMacros.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00010 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
00011 #include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"
00012 
00013 #include "FWCore/Framework/interface/ESHandle.h"
00014 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00015 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00016 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00017 
00018 
00019 PrimaryVertexProducer::PrimaryVertexProducer(const edm::ParameterSet& conf)
00020   :theConfig(conf)
00021 {
00022 
00023   fVerbose   = conf.getUntrackedParameter<bool>("verbose", false);
00024   trackLabel = conf.getParameter<edm::InputTag>("TrackLabel");
00025   beamSpotLabel = conf.getParameter<edm::InputTag>("beamSpotLabel");
00026 
00027 
00028   // select and configure the track selection
00029   std::string trackSelectionAlgorithm=conf.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
00030   if(trackSelectionAlgorithm=="filter"){
00031     theTrackFilter= new TrackFilterForPVFinding( conf.getParameter<edm::ParameterSet>("TkFilterParameters") );
00032   }else if (trackSelectionAlgorithm=="filterWithThreshold"){
00033     theTrackFilter= new HITrackFilterForPVFinding(conf.getParameter<edm::ParameterSet>("TkFilterParameters"));
00034   }else{
00035     throw VertexException("PrimaryVertexProducerAlgorithm: unknown track selection algorithm: " + trackSelectionAlgorithm);  
00036   }
00037 
00038 
00039   // select and configure the track clusterizer
00040   std::string clusteringAlgorithm=conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<std::string>("algorithm");
00041   if (clusteringAlgorithm=="gap"){
00042     theTrackClusterizer = new GapClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkGapClusParameters"));
00043   }else if(clusteringAlgorithm=="DA"){
00044     theTrackClusterizer = new DAClusterizerInZ(conf.getParameter<edm::ParameterSet>("TkClusParameters").getParameter<edm::ParameterSet>("TkDAClusParameters"));
00045   }else{
00046     throw VertexException("PrimaryVertexProducerAlgorithm: unknown clustering algorithm: " + clusteringAlgorithm);  
00047   }
00048 
00049 
00050   // select and configure the vertex fitters
00051   if (conf.exists("vertexCollections")){
00052     std::vector<edm::ParameterSet> vertexCollections =conf.getParameter< std::vector<edm::ParameterSet> >("vertexCollections");
00053 
00054     for( std::vector< edm::ParameterSet >::const_iterator algoconf = vertexCollections.begin(); algoconf != vertexCollections.end(); algoconf++){
00055       
00056       algo algorithm;
00057       std::string fitterAlgorithm = algoconf->getParameter<std::string>("algorithm");
00058       if (fitterAlgorithm=="KalmanVertexFitter") {
00059         algorithm.fitter= new KalmanVertexFitter();
00060       } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
00061         algorithm.fitter= new AdaptiveVertexFitter();
00062       } else {
00063         throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);  
00064       }
00065       algorithm.label = algoconf->getParameter<std::string>("label");
00066       algorithm.minNdof = algoconf->getParameter<double>("minNdof");
00067       algorithm.useBeamConstraint=algoconf->getParameter<bool>("useBeamConstraint");
00068       algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), algoconf->getParameter<double>("maxDistanceToBeam"));
00069       algorithms.push_back(algorithm);
00070       
00071       produces<reco::VertexCollection>(algorithm.label);
00072     }
00073   }else{
00074     edm::LogWarning("MisConfiguration")<<"this module's configuration has changed, please update to have a vertexCollections=cms.VPSet parameter.";
00075 
00076     algo algorithm;
00077     std::string fitterAlgorithm = conf.getParameter<std::string>("algorithm");
00078     if (fitterAlgorithm=="KalmanVertexFitter") {
00079       algorithm.fitter= new KalmanVertexFitter();
00080     } else if( fitterAlgorithm=="AdaptiveVertexFitter") {
00081       algorithm.fitter= new AdaptiveVertexFitter();
00082     } else {
00083       throw VertexException("PrimaryVertexProducerAlgorithm: unknown algorithm: " + fitterAlgorithm);  
00084     }
00085     algorithm.label = "";
00086     algorithm.minNdof = conf.getParameter<double>("minNdof");
00087     algorithm.useBeamConstraint=conf.getParameter<bool>("useBeamConstraint");
00088     
00089     algorithm.vertexSelector=new VertexCompatibleWithBeam(VertexDistanceXY(), conf.getParameter<edm::ParameterSet>("PVSelParameters").getParameter<double>("maxDistanceToBeam"));
00090 
00091     algorithms.push_back(algorithm);
00092     produces<reco::VertexCollection>(algorithm.label);
00093   }
00094  
00095 
00096 }
00097 
00098 
00099 PrimaryVertexProducer::~PrimaryVertexProducer() {
00100   if (theTrackFilter) delete theTrackFilter;
00101   if (theTrackClusterizer) delete theTrackClusterizer;
00102   for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
00103     if (algorithm->fitter) delete algorithm->fitter;
00104     if (algorithm->vertexSelector) delete algorithm->vertexSelector;
00105   }
00106 }
00107 
00108 
00109 void
00110 PrimaryVertexProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00111 {
00112 
00113   // get the BeamSpot, it will alwys be needed, even when not used as a constraint
00114   reco::BeamSpot beamSpot;
00115   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00116   iEvent.getByLabel(beamSpotLabel,recoBeamSpotHandle);
00117   if (recoBeamSpotHandle.isValid()){
00118     beamSpot = *recoBeamSpotHandle;
00119   }else{
00120     edm::LogError("UnusableBeamSpot") << "No beam spot available from EventSetup";
00121   }
00122 
00123   bool validBS = true;
00124   VertexState beamVertexState(beamSpot);
00125   if ( (beamVertexState.error().cxx() <= 0.) || 
00126        (beamVertexState.error().cyy() <= 0.) ||
00127        (beamVertexState.error().czz() <= 0.) ) {
00128     validBS = false;
00129     edm::LogError("UnusableBeamSpot") << "Beamspot with invalid errors "<<beamVertexState.error().matrix();
00130   }
00131 
00132 
00133   // get RECO tracks from the event
00134   // `tks` can be used as a ptr to a reco::TrackCollection
00135   edm::Handle<reco::TrackCollection> tks;
00136   iEvent.getByLabel(trackLabel, tks);
00137 
00138 
00139   // interface RECO tracks to vertex reconstruction
00140   edm::ESHandle<TransientTrackBuilder> theB;
00141   iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theB);
00142   std::vector<reco::TransientTrack> t_tks = (*theB).build(tks, beamSpot);
00143   if(fVerbose) {std::cout << "RecoVertex/PrimaryVertexProducer"
00144                      << "Found: " << t_tks.size() << " reconstructed tracks" << "\n";
00145   }
00146 
00147 
00148   // select tracks
00149   std::vector<reco::TransientTrack> seltks = theTrackFilter->select( t_tks );
00150 
00151 
00152   // clusterize tracks in Z
00153   std::vector< std::vector<reco::TransientTrack> > clusters =  theTrackClusterizer->clusterize(seltks);
00154   if (fVerbose){std::cout <<  " clustering returned  "<< clusters.size() << " clusters  from " << seltks.size() << " selected tracks" <<std::endl;}
00155 
00156 
00157   // vertex fits
00158   for( std::vector <algo>::const_iterator algorithm=algorithms.begin(); algorithm!=algorithms.end(); algorithm++){
00159 
00160 
00161     std::auto_ptr<reco::VertexCollection> result(new reco::VertexCollection);
00162     reco::VertexCollection vColl;
00163 
00164 
00165     std::vector<TransientVertex> pvs;
00166     for (std::vector< std::vector<reco::TransientTrack> >::const_iterator iclus
00167            = clusters.begin(); iclus != clusters.end(); iclus++) {
00168 
00169 
00170       TransientVertex v; 
00171       if( algorithm->useBeamConstraint && validBS &&((*iclus).size()>1) ){
00172         
00173         v = algorithm->fitter->vertex(*iclus, beamSpot);
00174         
00175       }else if( !(algorithm->useBeamConstraint) && ((*iclus).size()>1) ) {
00176       
00177         v = algorithm->fitter->vertex(*iclus); 
00178         
00179       }// else: no fit ==> v.isValid()=False
00180 
00181 
00182       if (fVerbose){
00183         if (v.isValid()) std::cout << "x,y,z=" << v.position().x() <<" " << v.position().y() << " " <<  v.position().z() << std::endl;
00184         else std::cout <<"Invalid fitted vertex\n";
00185       }
00186 
00187       if (v.isValid() 
00188             && (v.degreesOfFreedom()>=algorithm->minNdof) 
00189           && (!validBS || (*(algorithm->vertexSelector))(v,beamVertexState))
00190           ) pvs.push_back(v);
00191     }// end of cluster loop
00192 
00193     if(fVerbose){
00194       std::cout << "PrimaryVertexProducerAlgorithm::vertices  candidates =" << pvs.size() << std::endl;
00195     }
00196 
00197 
00198 
00199     
00200 
00201     // sort vertices by pt**2  vertex (aka signal vertex tagging)
00202     if(pvs.size()>1){
00203       sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
00204     }
00205 
00206 
00207 
00208     // convert transient vertices returned by the theAlgo to (reco) vertices
00209     for (std::vector<TransientVertex>::const_iterator iv = pvs.begin();
00210          iv != pvs.end(); iv++) {
00211       reco::Vertex v = *iv;
00212       vColl.push_back(v);
00213     }
00214 
00215     if (vColl.empty()) {
00216       GlobalError bse(beamSpot.rotatedCovariance3D());
00217       if ( (bse.cxx() <= 0.) || 
00218            (bse.cyy() <= 0.) ||
00219            (bse.czz() <= 0.) ) {
00220         AlgebraicSymMatrix33 we;
00221         we(0,0)=10000; we(1,1)=10000; we(2,2)=10000;
00222         vColl.push_back(reco::Vertex(beamSpot.position(), we,0.,0.,0));
00223         if(fVerbose){
00224           std::cout <<"RecoVertex/PrimaryVertexProducer: "
00225                     << "Beamspot with invalid errors "<<bse.matrix()<<std::endl;
00226           std::cout << "Will put Vertex derived from dummy-fake BeamSpot into Event.\n";
00227         }
00228       } else {
00229         vColl.push_back(reco::Vertex(beamSpot.position(), 
00230                                      beamSpot.rotatedCovariance3D(),0.,0.,0));
00231         if(fVerbose){
00232           std::cout <<"RecoVertex/PrimaryVertexProducer: "
00233                     << " will put Vertex derived from BeamSpot into Event.\n";
00234         }
00235       }
00236     }
00237 
00238     if(fVerbose){
00239       int ivtx=0;
00240       for(reco::VertexCollection::const_iterator v=vColl.begin(); 
00241           v!=vColl.end(); ++v){
00242         std::cout << "recvtx "<< ivtx++ 
00243                   << "#trk " << std::setw(3) << v->tracksSize()
00244                   << " chi2 " << std::setw(4) << v->chi2() 
00245                   << " ndof " << std::setw(3) << v->ndof() 
00246                   << " x "  << std::setw(6) << v->position().x() 
00247                   << " dx " << std::setw(6) << v->xError()
00248                   << " y "  << std::setw(6) << v->position().y() 
00249                   << " dy " << std::setw(6) << v->yError()
00250                   << " z "  << std::setw(6) << v->position().z() 
00251                   << " dz " << std::setw(6) << v->zError()
00252                   << std::endl;
00253       }
00254     }
00255 
00256   
00257     *result = vColl;
00258     iEvent.put(result, algorithm->label); 
00259   }
00260   
00261 }
00262 
00263 
00264 //define this as a plug-in
00265 DEFINE_FWK_MODULE(PrimaryVertexProducer);