CMS 3D CMS Logo

GeneralTracksImporterWithVeto.cc
Go to the documentation of this file.
11 
13 public:
15  edm::ConsumesCollector& sumes) :
16  BlockElementImporterBase(conf,sumes),
17  src_(sumes.consumes<reco::PFRecTrackCollection>(conf.getParameter<edm::InputTag>("source"))),
18  veto_(sumes.consumes<reco::PFRecTrackCollection>(conf.getParameter<edm::InputTag>("veto"))),
19  muons_(sumes.consumes<reco::MuonCollection>(conf.getParameter<edm::InputTag>("muonSrc"))),
20  DPtovPtCut_(conf.getParameter<std::vector<double> >("DPtOverPtCuts_byTrackAlgo")),
21  NHitCut_(conf.getParameter<std::vector<unsigned> >("NHitCuts_byTrackAlgo")),
22  useIterTracking_(conf.getParameter<bool>("useIterativeTracking")),
23  cleanBadConvBrems_(conf.existsAs<bool>("cleanBadConvertedBrems") ? conf.getParameter<bool>("cleanBadConvertedBrems") : false),
24  debug_(conf.getUntrackedParameter<bool>("debug",false)) {
25 
26  pfmu_ = std::unique_ptr<PFMuonAlgo>(new PFMuonAlgo(conf));
27 
28  }
29 
30  void importToBlock( const edm::Event& ,
31  ElementList& ) const override;
32 
33 private:
34  int muAssocToTrack( const reco::TrackRef& trackref,
35  const edm::Handle<reco::MuonCollection>& muonh) const;
36 
39  const std::vector<double> DPtovPtCut_;
40  const std::vector<unsigned> NHitCut_;
42 
43  std::unique_ptr<PFMuonAlgo> pfmu_;
44 
45 };
46 
49  "GeneralTracksImporterWithVeto");
50 
55  auto tracks = e.getHandle(src_);
56  auto vetosH = e.getHandle(veto_);
57  const auto& vetos = *vetosH;
58  std::unordered_set<unsigned> vetoed;
59  for(unsigned i = 0; i < vetos.size(); ++i ) {
60  vetoed.insert(vetos[i].trackRef().key());
61  }
62  auto muons = e.getHandle(muons_);
63  elems.reserve(elems.size() + tracks->size());
64  std::vector<bool> mask(tracks->size(),true);
65  reco::MuonRef muonref;
66 
67  // remove converted brems with bad pT resolution if requested
68  // this reproduces the old behavior of PFBlockAlgo
69  if( cleanBadConvBrems_ ) {
70  auto itr = elems.begin();
71  while( itr != elems.end() ) {
72  if( (*itr)->type() == reco::PFBlockElement::TRACK ) {
73  const reco::PFBlockElementTrack* trkel =
74  static_cast<reco::PFBlockElementTrack*>(itr->get());
75  const reco::ConversionRefVector& cRef = trkel->convRefs();
76  const reco::PFDisplacedTrackerVertexRef& dvRef =
79  trkel->V0Ref();
80  // if there is no displaced vertex reference and it is marked
81  // as a conversion it's gotta be a converted brem
83  cRef.empty() && dvRef.isNull() && v0Ref.isNull() ) {
84  // if the Pt resolution is bad we kill this element
86  itr = elems.erase(itr);
87  continue;
88  }
89  }
90  }
91  ++itr;
92  } // loop on existing elements
93  }
94  // preprocess existing tracks in the element list and create a mask
95  // so that we do not import tracks twice, tag muons we find
96  // in this collection
97  auto TKs_end = std::partition(elems.begin(),elems.end(),
98  [](const ElementType& a){
99  return a->type() == reco::PFBlockElement::TRACK;
100  });
101  auto btk_elems = elems.begin();
102  auto btrack = tracks->cbegin();
103  auto etrack = tracks->cend();
104  for( auto track = btrack; track != etrack; ++track) {
105  auto tk_elem = std::find_if(btk_elems,TKs_end,
106  [&](const ElementType& a){
107  return ( a->trackRef() ==
108  track->trackRef() );
109  });
110  if( tk_elem != TKs_end ) {
111  mask[std::distance(tracks->cbegin(),track)] = false;
112  // check and update if this track is a muon
113  const int muId = muAssocToTrack( (*tk_elem)->trackRef(), muons );
114  if( muId != -1 ) {
115  muonref= reco::MuonRef( muons, muId );
116  if( PFMuonAlgo::isLooseMuon(muonref) || PFMuonAlgo::isMuon(muonref) ) {
117  static_cast<reco::PFBlockElementTrack*>(tk_elem->get())->setMuonRef(muonref);
118  }
119  }
120  }
121  }
122  // now we actually insert tracks, again tagging muons along the way
123  reco::PFRecTrackRef pftrackref;
124  reco::PFBlockElementTrack* trkElem = nullptr;
125  for( auto track = btrack; track != etrack; ++track) {
126  const unsigned idx = std::distance(btrack,track);
127  // since we already set muon refs in the previously imported tracks,
128  // here we can skip everything that is already imported
129  if( !mask[idx] ) continue;
130  muonref = reco::MuonRef();
131  pftrackref = reco::PFRecTrackRef(tracks,idx);
132  // Get the eventual muon associated to this track
133  const int muId = muAssocToTrack( pftrackref->trackRef(), muons );
134  bool thisIsAPotentialMuon = false;
135  if( muId != -1 ) {
136  muonref= reco::MuonRef( muons, muId );
137  thisIsAPotentialMuon = ( (pfmu_->hasValidTrack(muonref,true)&&PFMuonAlgo::isLooseMuon(muonref)) ||
138  (pfmu_->hasValidTrack(muonref,false)&&PFMuonAlgo::isMuon(muonref)));
139  }
140  if(thisIsAPotentialMuon || PFTrackAlgoTools::goodPtResolution( pftrackref->trackRef(), DPtovPtCut_, NHitCut_, useIterTracking_, debug_ ) ) {
141  trkElem = new reco::PFBlockElementTrack( pftrackref );
142  if (thisIsAPotentialMuon && debug_) {
143  std::cout << "Potential Muon P " << pftrackref->trackRef()->p()
144  << " pt " << pftrackref->trackRef()->p() << std::endl;
145  }
146  if( muId != -1 ) trkElem->setMuonRef(muonref);
147  if( vetoed.count(pftrackref->trackRef().key()) == 0 || muonref.isNonnull()){
148  elems.emplace_back(trkElem);
149  }
150  else delete trkElem;
151  }
152  }
153  elems.shrink_to_fit();
154 }
155 
158  const edm::Handle<reco::MuonCollection>& muonh) const {
159  auto muon = std::find_if(muonh->cbegin(),muonh->cend(),
160  [&](const reco::Muon& m) {
161  return ( m.track().isNonnull() &&
162  m.track() == trackref );
163  });
164  return ( muon != muonh->cend() ? std::distance(muonh->cbegin(),muon) : -1 );
165 }
int muAssocToTrack(const reco::TrackRef &trackref, const edm::Handle< reco::MuonCollection > &muonh) const
static bool isMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:66
edm::EDGetTokenT< reco::MuonCollection > muons_
const reco::TrackRef & trackRef() const override
const ConversionRefVector & convRefs() const override
bool trackType(TrackType trType) const override
edm::EDGetTokenT< reco::PFRecTrackCollection > veto_
key_type key() const
Accessor for product key.
Definition: Ref.h:263
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
Handle< PROD > getHandle(EDGetTokenT< PROD > token) const
Definition: Event.h:539
edm::EDGetTokenT< reco::PFRecTrackCollection > src_
edm::Ref< PFRecTrackCollection > PFRecTrackRef
persistent reference to PFRecTrack objects
Definition: PFRecTrackFwd.h:15
const VertexCompositeCandidateRef & V0Ref() const override
edm::Ref< MuonCollection > MuonRef
presistent reference to a Muon
Definition: MuonFwd.h:13
bool isNull() const
Checks for null.
Definition: Ref.h:248
const PFDisplacedTrackerVertexRef & displacedVertexRef(TrackType trType) const override
GeneralTracksImporterWithVeto(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
void importToBlock(const edm::Event &, ElementList &) const override
static bool isLooseMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:78
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
void setMuonRef(const MuonRef &muref) override
reference to the Muon
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
bool goodPtResolution(const reco::TrackRef &, const std::vector< double > &DPtovPtCut, const std::vector< unsigned > &NHitCut, bool useIterTracking, bool debug=false)
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
const std::vector< unsigned > NHitCut_