Go to the documentation of this file.00001 #include "RecoTauTag/RecoTau/interface/RecoTauVertexAssociator.h"
00002
00003 #include <functional>
00004 #include <boost/foreach.hpp>
00005
00006 #include "DataFormats/TauReco/interface/PFTau.h"
00007 #include "DataFormats/VertexReco/interface/Vertex.h"
00008 #include "FWCore/Framework/interface/Event.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "RecoTauTag/RecoTau/interface/RecoTauCommonUtilities.h"
00011 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00012 #include "DataFormats/TrackReco/interface/Track.h"
00013
00014
00015 namespace reco { namespace tau {
00016
00017 bool vxTrkFiltering;
00018
00019
00020
00021 reco::TrackBaseRef RecoTauVertexAssociator::getLeadTrack(const PFJet& jet) const{
00022 std::vector<PFCandidatePtr> allTracks = pfChargedCands(jet, true);
00023 std::vector<PFCandidatePtr> tracks;
00024
00025 if (!allTracks.size()){
00026 LogDebug("VxTrkAssocInfo") << " No tracks at this jet! Returning empty reference.";
00027 return reco::TrackBaseRef();
00028 }
00029 if(vxTrkFiltering) tracks = qcuts_.filterRefs(allTracks);
00030 else{
00031 tracks = allTracks;
00032 LogDebug("VxTrkAssocInfo") << " No quality cuts applied. All tracks passing.";
00033 }
00034 PFCandidatePtr cand;
00035 if (!tracks.size()){
00036 if(!recoverLeadingTrk){
00037 LogDebug("VxTrkAssocInfo") << " All " << allTracks.size() << " tracks rejected!";
00038 return reco::TrackBaseRef();
00039 }else{
00040 cand = allTracks[0];
00041 LogDebug("VxTrkAssocInfo") << " All " << allTracks.size() << " tracks rejected but leading track was recovered!";
00042 }
00043 }else cand = tracks[0];
00044
00045 if (cand->trackRef().isNonnull())
00046 return reco::TrackBaseRef(cand->trackRef());
00047 else if (cand->gsfTrackRef().isNonnull()) {
00048 return reco::TrackBaseRef(cand->gsfTrackRef());
00049 }
00050 return reco::TrackBaseRef();
00051 }
00052 namespace {
00053
00054
00055 class DZtoTrack : public std::unary_function<double, reco::VertexRef> {
00056 public:
00057 DZtoTrack(const reco::TrackBaseRef& trk):trk_(trk){}
00058 double operator()(const reco::VertexRef& vtx) const {
00059 if (!trk_ || !vtx) {
00060 return std::numeric_limits<double>::infinity();
00061 }
00062 return std::abs(trk_->dz(vtx->position()));
00063 }
00064 private:
00065 const reco::TrackBaseRef trk_;
00066 };
00067
00068 class TrackWeightInVertex : public std::unary_function<double, reco::VertexRef>
00069 {
00070 public:
00071 TrackWeightInVertex(const reco::TrackBaseRef& trk):trk_(trk){}
00072 double operator()(const reco::VertexRef& vtx) const {
00073 if (!trk_ || !vtx) {
00074 return 0.0;
00075 }
00076 return vtx->trackWeight(trk_);
00077 }
00078 private:
00079 const reco::TrackBaseRef trk_;
00080 };
00081
00082 }
00083
00084 RecoTauVertexAssociator::RecoTauVertexAssociator(
00085 const edm::ParameterSet& pset): qcuts_(pset.exists("vxAssocQualityCuts") ? pset.getParameterSet("vxAssocQualityCuts") : pset.getParameterSet("signalQualityCuts"))
00086 {
00087
00088 vertexTag_ = edm::InputTag("offlinePrimaryVertices", "");
00089 std::string algorithm = "highestPtInEvent";
00090
00091 if (!pset.exists("primaryVertexSrc") || !pset.exists("pvFindingAlgo")) {
00092 edm::LogWarning("NoVertexFindingMethodSpecified")
00093 << "The PSet passed to the RecoTauVertexAssociator was"
00094 << " incorrectly configured. The vertex will be taken as the "
00095 << "highest Pt vertex from the offlinePrimaryVertices collection."
00096 << std::endl;
00097 } else {
00098 vertexTag_ = pset.getParameter<edm::InputTag>("primaryVertexSrc");
00099 algorithm = pset.getParameter<std::string>("pvFindingAlgo");
00100 }
00101 vxTrkFiltering = false;
00102 if(!pset.exists("vertexTrackFiltering") && pset.exists("vxAssocQualityCuts")){
00103 edm::LogWarning("NoVertexTrackFilteringSpecified")
00104 << "The PSet passed to the RecoTauVertexAssociator was"
00105 << " incorrectly configured. Please define vertexTrackFiltering in config file."
00106 << " No filtering of tracks to vertices will be applied"
00107 << std::endl;
00108 }else{
00109 vxTrkFiltering = pset.exists("vertexTrackFiltering") ? pset.getParameter<bool>("vertexTrackFiltering") : false;
00110 }
00111
00112 LogDebug("TauVxAssociatorInfo") << " ----> Using " << algorithm << " algorithm to select vertex.";
00113 if (algorithm == "highestPtInEvent") {
00114 algo_ = kHighestPtInEvent;
00115 } else if (algorithm == "closestInDeltaZ") {
00116 algo_ = kClosestDeltaZ;
00117 } else if (algorithm == "highestWeightForLeadTrack") {
00118 algo_ = kHighestWeigtForLeadTrack;
00119 } else if (algorithm == "combined"){
00120 algo_ = kCombined;
00121 } else {
00122 throw cms::Exception("BadVertexAssociatorConfig")
00123 << "The algorithm specified for tau-vertex association "
00124 << algorithm << " is invalid. Options are: " << std::endl
00125 << "highestPtInEvent, "
00126 << "closestInDeltaZ, "
00127 << "highestWeightForLeadTrack, "
00128 << " or combined." << std::endl;
00129 }
00130 recoverLeadingTrk = pset.exists("recoverLeadingTrk") ? pset.getParameter<bool>("recoverLeadingTrk") : false;
00131
00132 JetToVertexAssociation=0;
00133 myEventNumber = -999;
00134
00135 }
00136
00137
00138 void RecoTauVertexAssociator::setEvent(const edm::Event& evt) {
00139 edm::Handle<reco::VertexCollection> verticesH_;
00140 evt.getByLabel(vertexTag_, verticesH_);
00141 vertices_.clear();
00142 vertices_.reserve(verticesH_->size());
00143 for(size_t i = 0; i < verticesH_->size(); ++i) {
00144 vertices_.push_back(reco::VertexRef(verticesH_, i));
00145 }
00146 if(vertices_.size()>0 ) qcuts_.setPV(vertices_[0]);
00147 int currentEvent = evt.id().event();
00148 if(myEventNumber == -999 || myEventNumber!=currentEvent)
00149 {
00150 if(myEventNumber==-999) JetToVertexAssociation= new std::map<const reco::PFJet*,reco::VertexRef>;
00151 else JetToVertexAssociation->clear();
00152 myEventNumber = currentEvent;
00153 }
00154 }
00155
00156 reco::VertexRef
00157 RecoTauVertexAssociator::associatedVertex(const PFTau& tau) const {
00158 reco::PFJetRef jetRef = tau.jetRef();
00159
00160 if (jetRef.isNull())
00161 jetRef = tau.pfTauTagInfoRef()->pfjetRef();
00162 return associatedVertex(*jetRef);
00163 }
00164
00165 reco::VertexRef
00166 RecoTauVertexAssociator::associatedVertex(const PFJet& jet) const {
00167 reco::VertexRef output = vertices_.size() ? vertices_[0] : reco::VertexRef();
00168 PFJet const* my_jet_ptr = &jet;
00169 LogDebug("VxTrkAssocInfo") << "There are " << vertices_.size() << " vertices in this event. The jet is " << jet;
00170 LogTrace("VxTrkAssocInfo") << "The lenght of assoc map is "<< JetToVertexAssociation->size();
00171
00172
00173 std::map<const reco::PFJet*,reco::VertexRef>::iterator it = JetToVertexAssociation->find(my_jet_ptr);
00174 if(it!=JetToVertexAssociation->end())
00175 {
00176 LogTrace("VxTrkAssocInfo") << "I have seen this jet! Returning pointer to stored value.";
00177 return it->second;
00178 }
00179
00180 reco::TrackBaseRef leadTrack;
00181 if (algo_ == kHighestPtInEvent) {
00182 return output;
00183 } else if (algo_ == kClosestDeltaZ) {
00184 leadTrack = getLeadTrack(jet);
00185 if(!leadTrack){
00186 LogTrace("VxTrkAssocInfo") << "No track to associate to! Returning vx no. 0.";
00187 JetToVertexAssociation->insert(std::pair<const PFJet*, reco::VertexRef>(my_jet_ptr,output));
00188 return output;
00189 }
00190 double closestDistance = std::numeric_limits<double>::infinity();
00191 DZtoTrack dzComputer(leadTrack);
00192
00193 BOOST_FOREACH(const reco::VertexRef& vtx, vertices_) {
00194 double dz = dzComputer(vtx);
00195 if (dz < closestDistance) {
00196 closestDistance = dz;
00197 output = vtx;
00198
00199 }
00200 }
00201 } else if (algo_ == kHighestWeigtForLeadTrack || algo_ == kCombined) {
00202 leadTrack = getLeadTrack(jet);
00203 if(!leadTrack){
00204 LogTrace("VxTrkAssocInfo") << "No track to associate to! Returning vx no. 0.";
00205 JetToVertexAssociation->insert(std::pair<const PFJet*, reco::VertexRef>(my_jet_ptr,output));
00206 return output;
00207 }
00208 double largestWeight = 0.;
00209
00210 TrackWeightInVertex weightComputer(leadTrack);
00211 BOOST_FOREACH(const reco::VertexRef& vtx, vertices_) {
00212 double weight = weightComputer(vtx);
00213 if (weight > largestWeight) {
00214 largestWeight = weight;
00215 output = vtx;
00216 }
00217 }
00218
00219 if(algo_==kCombined && largestWeight < 1e-7){
00220 LogTrace("VxTrkAssocInfo") << " No vertex had positive weight! Trying dZ instead... ";
00221 DZtoTrack dzComputer(leadTrack);
00222 double closestDistance = std::numeric_limits<double>::infinity();
00223 BOOST_FOREACH(const reco::VertexRef& vtx, vertices_) {
00224 double dz_i = dzComputer(vtx);
00225 if (dz_i < closestDistance) {
00226 closestDistance = dz_i;
00227 output = vtx;
00228 }
00229 }
00230 }
00231 }
00232
00233 JetToVertexAssociation->insert(std::pair<const PFJet*, reco::VertexRef>(my_jet_ptr,output));
00234
00235 return output;
00236 }
00237
00238
00239
00240 }}