CMS 3D CMS Logo

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