Go to the documentation of this file.00001 #include <memory>
00002 #include <set>
00003
00004
00005 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00006 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00007 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00008 #include "TrackingTools/IPTools/interface/IPTools.h"
00009
00010
00011 #include "FWCore/Framework/interface/EDProducer.h"
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/MakerMacros.h"
00014 #include "FWCore/Utilities/interface/InputTag.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016
00017 #include "DataFormats/Common/interface/Handle.h"
00018 #include "DataFormats/TrackReco/interface/Track.h"
00019 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00020 #include "DataFormats/VertexReco/interface/Vertex.h"
00021 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00022 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
00023
00024
00025 #include "DataFormats/TrackReco/interface/Track.h"
00026 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00027 #include "DataFormats/VertexReco/interface/Vertex.h"
00028 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00029 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00030
00031 #include "RecoVertex/AdaptiveVertexFit/interface/AdaptiveVertexFitter.h"
00032 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexUpdator.h"
00033 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexTrackCompatibilityEstimator.h"
00034 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexSmoother.h"
00035 #include "DataFormats/Math/interface/deltaR.h"
00036
00037
00038 #include "RecoVertex/ConfigurableVertexReco/interface/ConfigurableVertexReconstructor.h"
00039
00040
00041
00042 class TrackVertexArbitrator : public edm::EDProducer {
00043 public:
00044 TrackVertexArbitrator(const edm::ParameterSet ¶ms);
00045
00046
00047 virtual void produce(edm::Event &event, const edm::EventSetup &es);
00048
00049 private:
00050 bool trackFilter(const reco::TrackRef &track) const;
00051
00052 edm::InputTag primaryVertexCollection;
00053 edm::InputTag secondaryVertexCollection;
00054 edm::InputTag trackCollection;
00055 edm::InputTag beamSpotCollection;
00056 double dRCut;
00057 double distCut;
00058 double sigCut;
00059 double dLenFraction;
00060 };
00061
00062 TrackVertexArbitrator::TrackVertexArbitrator(const edm::ParameterSet ¶ms) :
00063 primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")),
00064 secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")),
00065 trackCollection(params.getParameter<edm::InputTag>("tracks")),
00066 beamSpotCollection(params.getParameter<edm::InputTag>("beamSpot")),
00067 dRCut(params.getParameter<double>("dRCut")),
00068 distCut(params.getParameter<double>("distCut")),
00069 sigCut(params.getParameter<double>("sigCut")),
00070 dLenFraction(params.getParameter<double>("dLenFraction"))
00071 {
00072 produces<reco::VertexCollection>();
00073 }
00074
00075 bool TrackVertexArbitrator::trackFilter(const reco::TrackRef &track) const
00076 {
00077 if (track->hitPattern().trackerLayersWithMeasurement() < 4)
00078 return false;
00079 if (track->pt() < 0.4 )
00080 return false;
00081
00082 return true;
00083 }
00084
00085
00086 void TrackVertexArbitrator::produce(edm::Event &event, const edm::EventSetup &es)
00087 {
00088 using namespace reco;
00089
00090 edm::Handle<VertexCollection> secondaryVertices;
00091 event.getByLabel(secondaryVertexCollection, secondaryVertices);
00092
00093 edm::Handle<VertexCollection> primaryVertices;
00094 event.getByLabel(primaryVertexCollection, primaryVertices);
00095
00096 edm::Handle<TrackCollection> tracks;
00097 event.getByLabel(trackCollection, tracks);
00098
00099 edm::ESHandle<TransientTrackBuilder> trackBuilder;
00100 es.get<TransientTrackRecord>().get("TransientTrackBuilder",
00101 trackBuilder);
00102
00103 edm::Handle<BeamSpot> beamSpot;
00104 event.getByLabel(beamSpotCollection, beamSpot);
00105
00106 const reco::Vertex &pv = (*primaryVertices)[0];
00107
00108 VertexDistance3D dist;
00109
00110 double sigmacut = 3.0;
00111 double Tini = 256.;
00112 double ratio = 0.25;
00113
00114 AdaptiveVertexFitter theAdaptiveFitter(
00115 GeometricAnnealing(sigmacut, Tini, ratio),
00116 DefaultLinearizationPointFinder(),
00117 KalmanVertexUpdator<5>(),
00118 KalmanVertexTrackCompatibilityEstimator<5>(),
00119 KalmanVertexSmoother() );
00120
00121
00122
00123 std::auto_ptr<VertexCollection> recoVertices(new VertexCollection);
00124
00125 VertexDistance3D vdist;
00126
00127 for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin();
00128 sv != secondaryVertices->end(); ++sv) {
00129
00130
00131
00132
00133
00134
00135 GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z());
00136 GlobalPoint ssv(sv->position().x(),sv->position().y(),sv->position().z());
00137 GlobalVector flightDir = ssv-ppv;
00138
00139 Measurement1D dlen= vdist.distance(pv,*sv);
00140 std::vector<reco::TransientTrack> selTracks;
00141
00142 for(TrackCollection::const_iterator track = tracks->begin();
00143 track != tracks->end(); ++track) {
00144
00145 TrackRef ref(tracks, track - tracks->begin());
00146 if (!trackFilter(ref)) continue;
00147
00148 TransientTrack tt = trackBuilder->build(ref);
00149 tt.setBeamSpot(*beamSpot);
00150 float w = sv->trackWeight(ref);
00151 std::pair<bool,Measurement1D> ipv = IPTools::absoluteImpactParameter3D(tt,pv);
00152 std::pair<bool,Measurement1D> isv = IPTools::absoluteImpactParameter3D(tt,*sv);
00153 if( w > 0 || ( isv.second.significance() < sigCut && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction ) )
00154 {
00155 float dR = deltaR(flightDir.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi());
00156
00157 if(isv.second.value() < ipv.second.value() && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction
00158 && dR < dRCut )
00159 {
00160
00161
00162 selTracks.push_back(tt);
00163 } else
00164 {
00165
00166
00167
00168 if(w > 0.5 && isv.second.value() < ipv.second.value() && dR < dRCut) selTracks.push_back(tt);
00169 }
00170
00171
00172
00173
00174
00175 }
00176 }
00177
00178 if(selTracks.size() >= 2)
00179 {
00180 TransientVertex singleFitVertex;
00181 singleFitVertex = theAdaptiveFitter.vertex(selTracks,ssv);
00182 if(singleFitVertex.isValid()) recoVertices->push_back(singleFitVertex);
00183 }
00184 }
00185 event.put(recoVertices);
00186 }
00187
00188 DEFINE_FWK_MODULE(TrackVertexArbitrator);