Public Member Functions | |
virtual void | produce (edm::Event &event, const edm::EventSetup &es) |
TrackVertexArbitrator (const edm::ParameterSet ¶ms) | |
Private Member Functions | |
bool | trackFilter (const reco::TrackRef &track) const |
Private Attributes | |
edm::InputTag | beamSpotCollection |
double | distCut |
double | dLenFraction |
double | dRCut |
edm::InputTag | primaryVertexCollection |
edm::InputTag | secondaryVertexCollection |
double | sigCut |
edm::InputTag | trackCollection |
Definition at line 42 of file TrackVertexArbitrator.cc.
TrackVertexArbitrator::TrackVertexArbitrator | ( | const edm::ParameterSet & | params | ) |
Definition at line 62 of file TrackVertexArbitrator.cc.
: primaryVertexCollection(params.getParameter<edm::InputTag>("primaryVertices")), secondaryVertexCollection(params.getParameter<edm::InputTag>("secondaryVertices")), trackCollection(params.getParameter<edm::InputTag>("tracks")), beamSpotCollection(params.getParameter<edm::InputTag>("beamSpot")), dRCut(params.getParameter<double>("dRCut")), distCut(params.getParameter<double>("distCut")), sigCut(params.getParameter<double>("sigCut")), dLenFraction(params.getParameter<double>("dLenFraction")) { produces<reco::VertexCollection>(); }
void TrackVertexArbitrator::produce | ( | edm::Event & | event, |
const edm::EventSetup & | es | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 86 of file TrackVertexArbitrator.cc.
References IPTools::absoluteImpactParameter3D(), SiPixelRawToDigiRegional_cfi::beamSpot, beamSpotCollection, deltaR(), distCut, dLenFraction, dRCut, reco::TrackBase::eta(), edm::EventSetup::get(), TransientVertex::isValid(), reco::TrackBase::phi(), primaryVertexCollection, dt_dqm_sourceclient_common_cff::reco, secondaryVertexCollection, reco::TransientTrack::setBeamSpot(), sigCut, reco::TransientTrack::track(), trackCollection, trackFilter(), testEve_cfg::tracks, groupFilesInBlocks::tt, Measurement1D::value(), AdaptiveVertexFitter::vertex(), and w().
{ using namespace reco; edm::Handle<VertexCollection> secondaryVertices; event.getByLabel(secondaryVertexCollection, secondaryVertices); edm::Handle<VertexCollection> primaryVertices; event.getByLabel(primaryVertexCollection, primaryVertices); edm::Handle<TrackCollection> tracks; event.getByLabel(trackCollection, tracks); edm::ESHandle<TransientTrackBuilder> trackBuilder; es.get<TransientTrackRecord>().get("TransientTrackBuilder", trackBuilder); edm::Handle<BeamSpot> beamSpot; event.getByLabel(beamSpotCollection, beamSpot); const reco::Vertex &pv = (*primaryVertices)[0]; // std::cout << "PV: " << pv.position() << std::endl; VertexDistance3D dist; double sigmacut = 3.0; double Tini = 256.; double ratio = 0.25; AdaptiveVertexFitter theAdaptiveFitter( GeometricAnnealing(sigmacut, Tini, ratio), DefaultLinearizationPointFinder(), KalmanVertexUpdator<5>(), KalmanVertexTrackCompatibilityEstimator<5>(), KalmanVertexSmoother() ); std::auto_ptr<VertexCollection> recoVertices(new VertexCollection); VertexDistance3D vdist; for(std::vector<reco::Vertex>::const_iterator sv = secondaryVertices->begin(); sv != secondaryVertices->end(); ++sv) { /* recoVertices->push_back(*sv); for(std::vector<reco::Vertex>::iterator sv = recoVertices->begin(); sv != recoVertices->end(); ++sv) { */ GlobalPoint ppv(pv.position().x(),pv.position().y(),pv.position().z()); GlobalPoint ssv(sv->position().x(),sv->position().y(),sv->position().z()); GlobalVector flightDir = ssv-ppv; // std::cout << "Vertex : " << sv-secondaryVertices->begin() << " " << sv->position() << std::endl; Measurement1D dlen= vdist.distance(pv,*sv); std::vector<reco::TransientTrack> selTracks; for(TrackCollection::const_iterator track = tracks->begin(); track != tracks->end(); ++track) { TrackRef ref(tracks, track - tracks->begin()); if (!trackFilter(ref)) continue; TransientTrack tt = trackBuilder->build(ref); tt.setBeamSpot(*beamSpot); float w = sv->trackWeight(ref); std::pair<bool,Measurement1D> ipv = IPTools::absoluteImpactParameter3D(tt,pv); std::pair<bool,Measurement1D> isv = IPTools::absoluteImpactParameter3D(tt,*sv); if( w > 0 || ( isv.second.significance() < sigCut && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction ) ) { float dR = deltaR(flightDir.eta(), flightDir.phi(), tt.track().eta(), tt.track().phi()); if(isv.second.value() < ipv.second.value() && isv.second.value() < distCut && isv.second.value() < dlen.value()*dLenFraction && dR < dRCut ) { // if(w > 0.5) std::cout << " = "; // else std::cout << " + "; selTracks.push_back(tt); } else { // if(w > 0.5 && isv.second.value() > ipv.second.value() ) std::cout << " - "; // else std::cout << " "; //add also the tracks used in previous fitting that are still closer to Sv than Pv if(w > 0.5 && isv.second.value() < ipv.second.value() && dR < dRCut) selTracks.push_back(tt); } // std::cout << "t : " << track-tracks->begin() << " w: " << w // << " svip: " << isv.second.significance() << " " << isv.second.value() // << " pvip: " << ipv.second.significance() << " " << ipv.second.value() << " dr: " << dR << std::endl; } } if(selTracks.size() >= 2) { TransientVertex singleFitVertex; singleFitVertex = theAdaptiveFitter.vertex(selTracks,ssv); if(singleFitVertex.isValid()) recoVertices->push_back(singleFitVertex); } } event.put(recoVertices); }
bool TrackVertexArbitrator::trackFilter | ( | const reco::TrackRef & | track | ) | const [private] |
Definition at line 75 of file TrackVertexArbitrator.cc.
Referenced by produce().
{ if (track->hitPattern().trackerLayersWithMeasurement() < 4) return false; if (track->pt() < 0.4 ) return false; return true; }
Definition at line 55 of file TrackVertexArbitrator.cc.
Referenced by produce().
double TrackVertexArbitrator::distCut [private] |
Definition at line 57 of file TrackVertexArbitrator.cc.
Referenced by produce().
double TrackVertexArbitrator::dLenFraction [private] |
Definition at line 59 of file TrackVertexArbitrator.cc.
Referenced by produce().
double TrackVertexArbitrator::dRCut [private] |
Definition at line 56 of file TrackVertexArbitrator.cc.
Referenced by produce().
Definition at line 52 of file TrackVertexArbitrator.cc.
Referenced by produce().
Definition at line 53 of file TrackVertexArbitrator.cc.
Referenced by produce().
double TrackVertexArbitrator::sigCut [private] |
Definition at line 58 of file TrackVertexArbitrator.cc.
Referenced by produce().
Definition at line 54 of file TrackVertexArbitrator.cc.
Referenced by produce().