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());
27  pfmu_->setParameters(conf);
28 
29  }
30 
31  void importToBlock( const edm::Event& ,
32  ElementList& ) const override;
33 
34 private:
35  int muAssocToTrack( const reco::TrackRef& trackref,
36  const edm::Handle<reco::MuonCollection>& muonh) const;
37 
40  const std::vector<double> DPtovPtCut_;
41  const std::vector<unsigned> NHitCut_;
43 
44  std::unique_ptr<PFMuonAlgo> pfmu_;
45 
46 };
47 
50  "GeneralTracksImporterWithVeto");
51 
57  e.getByToken(src_,tracks);
59  e.getByToken(veto_,vetosH);
60  const auto& vetos = *vetosH;
61  std::unordered_set<unsigned> vetoed;
62  for(unsigned i = 0; i < vetos.size(); ++i ) {
63  vetoed.insert(vetos[i].trackRef().key());
64  }
66  e.getByToken(muons_,muons);
67  elems.reserve(elems.size() + tracks->size());
68  std::vector<bool> mask(tracks->size(),true);
69  reco::MuonRef muonref;
70 
71  // remove converted brems with bad pT resolution if requested
72  // this reproduces the old behavior of PFBlockAlgo
73  if( cleanBadConvBrems_ ) {
74  auto itr = elems.begin();
75  while( itr != elems.end() ) {
76  if( (*itr)->type() == reco::PFBlockElement::TRACK ) {
77  const reco::PFBlockElementTrack* trkel =
78  static_cast<reco::PFBlockElementTrack*>(itr->get());
79  const reco::ConversionRefVector& cRef = trkel->convRefs();
80  const reco::PFDisplacedTrackerVertexRef& dvRef =
83  trkel->V0Ref();
84  // if there is no displaced vertex reference and it is marked
85  // as a conversion it's gotta be a converted brem
87  cRef.empty() && dvRef.isNull() && v0Ref.isNull() ) {
88  // if the Pt resolution is bad we kill this element
90  itr = elems.erase(itr);
91  continue;
92  }
93  }
94  }
95  ++itr;
96  } // loop on existing elements
97  }
98  // preprocess existing tracks in the element list and create a mask
99  // so that we do not import tracks twice, tag muons we find
100  // in this collection
101  auto TKs_end = std::partition(elems.begin(),elems.end(),
102  [](const ElementType& a){
103  return a->type() == reco::PFBlockElement::TRACK;
104  });
105  auto btk_elems = elems.begin();
106  auto btrack = tracks->cbegin();
107  auto etrack = tracks->cend();
108  for( auto track = btrack; track != etrack; ++track) {
109  auto tk_elem = std::find_if(btk_elems,TKs_end,
110  [&](const ElementType& a){
111  return ( a->trackRef() ==
112  track->trackRef() );
113  });
114  if( tk_elem != TKs_end ) {
115  mask[std::distance(tracks->cbegin(),track)] = false;
116  // check and update if this track is a muon
117  const int muId = muAssocToTrack( (*tk_elem)->trackRef(), muons );
118  if( muId != -1 ) {
119  muonref= reco::MuonRef( muons, muId );
120  if( PFMuonAlgo::isLooseMuon(muonref) || PFMuonAlgo::isMuon(muonref) ) {
121  static_cast<reco::PFBlockElementTrack*>(tk_elem->get())->setMuonRef(muonref);
122  }
123  }
124  }
125  }
126  // now we actually insert tracks, again tagging muons along the way
127  reco::PFRecTrackRef pftrackref;
128  reco::PFBlockElementTrack* trkElem = nullptr;
129  for( auto track = btrack; track != etrack; ++track) {
130  const unsigned idx = std::distance(btrack,track);
131  // since we already set muon refs in the previously imported tracks,
132  // here we can skip everything that is already imported
133  if( !mask[idx] ) continue;
134  muonref = reco::MuonRef();
135  pftrackref = reco::PFRecTrackRef(tracks,idx);
136  // Get the eventual muon associated to this track
137  const int muId = muAssocToTrack( pftrackref->trackRef(), muons );
138  bool thisIsAPotentialMuon = false;
139  if( muId != -1 ) {
140  muonref= reco::MuonRef( muons, muId );
141  thisIsAPotentialMuon = ( (pfmu_->hasValidTrack(muonref,true)&&PFMuonAlgo::isLooseMuon(muonref)) ||
142  (pfmu_->hasValidTrack(muonref,false)&&PFMuonAlgo::isMuon(muonref)));
143  }
144  if(thisIsAPotentialMuon || PFTrackAlgoTools::goodPtResolution( pftrackref->trackRef(), DPtovPtCut_, NHitCut_, useIterTracking_, debug_ ) ) {
145  trkElem = new reco::PFBlockElementTrack( pftrackref );
146  if (thisIsAPotentialMuon && debug_) {
147  std::cout << "Potential Muon P " << pftrackref->trackRef()->p()
148  << " pt " << pftrackref->trackRef()->p() << std::endl;
149  }
150  if( muId != -1 ) trkElem->setMuonRef(muonref);
151  if( vetoed.count(pftrackref->trackRef().key()) == 0 || muonref.isNonnull()){
152  elems.emplace_back(trkElem);
153  }
154  else delete trkElem;
155  }
156  }
157  elems.shrink_to_fit();
158 }
159 
162  const edm::Handle<reco::MuonCollection>& muonh) const {
163  auto muon = std::find_if(muonh->cbegin(),muonh->cend(),
164  [&](const reco::Muon& m) {
165  return ( m.track().isNonnull() &&
166  m.track() == trackref );
167  });
168  return ( muon != muonh->cend() ? std::distance(muonh->cbegin(),muon) : -1 );
169 }
int muAssocToTrack(const reco::TrackRef &trackref, const edm::Handle< reco::MuonCollection > &muonh) const
static bool isMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:155
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
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:265
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
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:250
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:168
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_