CMS 3D CMS Logo

GeneralTracksImporter.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  muons_(sumes.consumes<reco::MuonCollection>(conf.getParameter<edm::InputTag>("muonSrc"))),
19  DPtovPtCut_(conf.getParameter<std::vector<double> >("DPtOverPtCuts_byTrackAlgo")),
20  NHitCut_(conf.getParameter<std::vector<unsigned> >("NHitCuts_byTrackAlgo")),
21  useTiming_(conf.existsAs<edm::InputTag>("timeValueMap")),
22  srcTime_(useTiming_ ? sumes.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeValueMap")) : edm::EDGetTokenT<edm::ValueMap<float>>()),
23  srcTimeError_(useTiming_ ? sumes.consumes<edm::ValueMap<float>>(conf.getParameter<edm::InputTag>("timeErrorMap")) : edm::EDGetTokenT<edm::ValueMap<float>>()),
24  useIterTracking_(conf.getParameter<bool>("useIterativeTracking")),
25  cleanBadConvBrems_(conf.existsAs<bool>("cleanBadConvertedBrems") ? conf.getParameter<bool>("cleanBadConvertedBrems") : false),
26  debug_(conf.getUntrackedParameter<bool>("debug",false)) {
27 
28  pfmu_ = std::unique_ptr<PFMuonAlgo>(new PFMuonAlgo());
29  pfmu_->setParameters(conf);
30 
31  }
32 
33  void importToBlock( const edm::Event& ,
34  ElementList& ) const override;
35 
36 private:
37  int muAssocToTrack( const reco::TrackRef& trackref,
38  const edm::Handle<reco::MuonCollection>& muonh) const;
39 
42  const std::vector<double> DPtovPtCut_;
43  const std::vector<unsigned> NHitCut_;
44  const bool useTiming_;
47 
48  std::unique_ptr<PFMuonAlgo> pfmu_;
49 
50 };
51 
54  "GeneralTracksImporter");
55 
61  e.getByToken(src_,tracks);
63  e.getByToken(muons_,muons);
64  elems.reserve(elems.size() + tracks->size());
65  std::vector<bool> mask(tracks->size(),true);
66  reco::MuonRef muonref;
67  edm::Handle<edm::ValueMap<float>> timeH, timeErrH;
68  if (useTiming_) {
69  e.getByToken(srcTime_, timeH);
70  e.getByToken(srcTimeError_, timeErrH);
71  }
72  // remove converted brems with bad pT resolution if requested
73  // this reproduces the old behavior of PFBlockAlgo
74  if( cleanBadConvBrems_ ) {
75  auto itr = elems.begin();
76  while( itr != elems.end() ) {
77  if( (*itr)->type() == reco::PFBlockElement::TRACK ) {
78  const reco::PFBlockElementTrack* trkel =
79  static_cast<reco::PFBlockElementTrack*>(itr->get());
80  const reco::ConversionRefVector& cRef = trkel->convRefs();
81  const reco::PFDisplacedTrackerVertexRef& dvRef =
84  trkel->V0Ref();
85  // if there is no displaced vertex reference and it is marked
86  // as a conversion it's gotta be a converted brem
88  cRef.size() == 0 && dvRef.isNull() && v0Ref.isNull() ) {
89  // if the Pt resolution is bad we kill this element
91  itr = elems.erase(itr);
92  continue;
93  }
94  }
95  }
96  ++itr;
97  } // loop on existing elements
98  }
99  // preprocess existing tracks in the element list and create a mask
100  // so that we do not import tracks twice, tag muons we find
101  // in this collection
102  auto TKs_end = std::partition(elems.begin(),elems.end(),
103  [](const ElementType& a){
104  return a->type() == reco::PFBlockElement::TRACK;
105  });
106  auto btk_elems = elems.begin();
107  auto btrack = tracks->cbegin();
108  auto etrack = tracks->cend();
109  for( auto track = btrack; track != etrack; ++track) {
110  auto tk_elem = std::find_if(btk_elems,TKs_end,
111  [&](const ElementType& a){
112  return ( a->trackRef() ==
113  track->trackRef() );
114  });
115  if( tk_elem != TKs_end ) {
116  mask[std::distance(tracks->cbegin(),track)] = false;
117  // check and update if this track is a muon
118  const int muId = muAssocToTrack( (*tk_elem)->trackRef(), muons );
119  if( muId != -1 ) {
120  muonref= reco::MuonRef( muons, muId );
121  if( PFMuonAlgo::isLooseMuon(muonref) || PFMuonAlgo::isMuon(muonref) ) {
122  static_cast<reco::PFBlockElementTrack*>(tk_elem->get())->setMuonRef(muonref);
123  }
124  }
125  }
126  }
127  // now we actually insert tracks, again tagging muons along the way
128  reco::PFRecTrackRef pftrackref;
129  reco::PFBlockElementTrack* trkElem = NULL;
130  for( auto track = btrack; track != etrack; ++track) {
131  const unsigned idx = std::distance(btrack,track);
132  // since we already set muon refs in the previously imported tracks,
133  // here we can skip everything that is already imported
134  if( !mask[idx] ) continue;
135  muonref = reco::MuonRef();
136  pftrackref = reco::PFRecTrackRef(tracks,idx);
137  // Get the eventual muon associated to this track
138  const int muId = muAssocToTrack( pftrackref->trackRef(), muons );
139  bool thisIsAPotentialMuon = false;
140  if( muId != -1 ) {
141  muonref= reco::MuonRef( muons, muId );
142  thisIsAPotentialMuon = ( (pfmu_->hasValidTrack(muonref,true)&&PFMuonAlgo::isLooseMuon(muonref)) ||
143  (pfmu_->hasValidTrack(muonref,false)&&PFMuonAlgo::isMuon(muonref)));
144  }
145  if(thisIsAPotentialMuon || PFTrackAlgoTools::goodPtResolution( pftrackref->trackRef(), DPtovPtCut_, NHitCut_, useIterTracking_, debug_ ) ) {
146  trkElem = new reco::PFBlockElementTrack( pftrackref );
147  if (thisIsAPotentialMuon && debug_) {
148  std::cout << "Potential Muon P " << pftrackref->trackRef()->p()
149  << " pt " << pftrackref->trackRef()->p() << std::endl;
150  }
151  if( muId != -1 ) trkElem->setMuonRef(muonref);
152  if ( useTiming_ ) {
153  trkElem->setTime( (*timeH)[pftrackref->trackRef()], (*timeErrH)[pftrackref->trackRef()] );
154  }
155  elems.emplace_back(trkElem);
156  }
157  }
158  elems.shrink_to_fit();
159 }
160 
163  const edm::Handle<reco::MuonCollection>& muonh) const {
164  auto muon = std::find_if(muonh->cbegin(),muonh->cend(),
165  [&](const reco::Muon& m) {
166  return ( m.track().isNonnull() &&
167  m.track() == trackref );
168  });
169  return ( muon != muonh->cend() ? std::distance(muonh->cbegin(),muon) : -1 );
170 }
const VertexCompositeCandidateRef & V0Ref() const
void importToBlock(const edm::Event &, ElementList &) const override
const reco::TrackRef & trackRef() const
GeneralTracksImporter(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
edm::EDGetTokenT< edm::ValueMap< float > > srcTime_
const std::vector< double > DPtovPtCut_
edm::EDGetTokenT< reco::MuonCollection > muons_
static bool isMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:155
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
std::unique_ptr< PFMuonAlgo > pfmu_
#define NULL
Definition: scimark2.h:8
int muAssocToTrack(const reco::TrackRef &trackref, const edm::Handle< reco::MuonCollection > &muonh) const
edm::EDGetTokenT< reco::PFRecTrackCollection > src_
std::vector< Muon > MuonCollection
collection of Muon objects
Definition: MuonFwd.h:9
const ConversionRefVector & convRefs() const
edm::Ref< PFRecTrackCollection > PFRecTrackRef
persistent reference to PFRecTrack objects
Definition: PFRecTrackFwd.h:15
edm::Ref< MuonCollection > MuonRef
presistent reference to a Muon
Definition: MuonFwd.h:13
const std::vector< unsigned > NHitCut_
bool isNull() const
Checks for null.
Definition: Ref.h:249
virtual bool trackType(TrackType trType) const
void setMuonRef(const MuonRef &muref)
reference to the Muon
static bool isLooseMuon(const reco::PFBlockElement &elt)
Definition: PFMuonAlgo.cc:168
fixed size matrix
HLT enums.
double a
Definition: hdecay.h:121
const PFDisplacedTrackerVertexRef & displacedVertexRef(TrackType trType) const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::unique_ptr< reco::PFBlockElement > > ElementList
edm::EDGetTokenT< edm::ValueMap< float > > srcTimeError_
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
void setTime(float time, float timeError=0.f)
the timing information