CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes
HIMultiTrackSelector Class Reference

#include <HIMultiTrackSelector.h>

Inheritance diagram for HIMultiTrackSelector:
edm::stream::EDProducer<>

Public Member Functions

 HIMultiTrackSelector ()
 constructor More...
 
 HIMultiTrackSelector (const edm::ParameterSet &cfg)
 
 ~HIMultiTrackSelector () override
 destructor More...
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Protected Types

typedef math::XYZPoint Point
 

Protected Member Functions

void beginStream (edm::StreamID) final
 
void ParseForestVars ()
 
void processMVA (edm::Event &evt, const edm::EventSetup &es, std::vector< float > &mvaVals_, const reco::VertexCollection &hVtx) const
 
void produce (edm::Event &evt, const edm::EventSetup &es) final
 process one event More...
 
virtual void run (edm::Event &evt, const edm::EventSetup &es) const
 
bool select (unsigned tsNum, const reco::BeamSpot &vertexBeamSpot, const TrackingRecHitCollection &recHits, const reco::Track &tk, const std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr, double mvaVal) const
 return class, or -1 if rejected More...
 
void selectVertices (unsigned int tsNum, const reco::VertexCollection &vtxs, std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr) const
 

Protected Attributes

std::vector< bool > applyAbsCutsIfNoPV_
 
std::vector< bool > applyAdaptedPVCuts_
 
bool applyPixelMergingCuts_
 
edm::EDGetTokenT< reco::BeamSpotbeamspot_
 
std::vector< double > chi2n_no1Dmod_par_
 
std::vector< double > chi2n_par_
 
std::vector< std::vector< double > > d0_par1_
 
std::vector< std::vector< double > > d0_par2_
 
std::string dbFileName_
 
std::vector< std::vector< double > > dz_par1_
 
std::vector< std::vector< double > > dz_par2_
 
GBRForestforest_
 
std::string forestLabel_
 
edm::ESGetToken< GBRForest, GBRWrapperRcdforestToken_
 
std::vector< std::string > forestVars_
 
edm::EDGetTokenT< TrackingRecHitCollectionhSrc_
 
std::vector< bool > keepAllTracks_
 
std::vector< double > max_d0_
 Impact parameter absolute cuts. More...
 
std::vector< double > max_d0NoPV_
 
std::vector< double > max_eta_
 
std::vector< int32_t > max_lostHitFraction_
 
std::vector< uint32_t > max_lostLayers_
 
std::vector< int32_t > max_minMissHitOutOrIn_
 
std::vector< double > max_relpterr_
 
std::vector< double > max_z0_
 
std::vector< double > max_z0NoPV_
 
std::vector< uint32_t > min_3Dlayers_
 
std::vector< double > min_eta_
 
std::vector< uint32_t > min_hits_bypass_
 
std::vector< uint32_t > min_layers_
 Cuts on numbers of layers with hits/3D hits/lost hits. More...
 
std::vector< double > min_MVA_
 
std::vector< uint32_t > min_nhits_
 
std::string mvaType_
 
std::vector< int > mvavars_indices
 
std::vector< std::string > name_
 
std::vector< double > nSigmaZ_
 
std::vector< std::vector< double > > pixel_pTMaxCut_
 
std::vector< std::vector< double > > pixel_pTMinCut_
 
std::vector< unsigned int > preFilter_
 
std::vector< reco::TrackBase::TrackQualityqualityToSet_
 
std::vector< std::vector< double > > res_par_
 
std::vector< bool > setQualityBit_
 do I have to set a quality bit? More...
 
edm::EDGetTokenT< reco::TrackCollectionsrc_
 source collection label More...
 
bool useAnyMVA_
 
bool useForestFromDB_
 
std::vector< bool > useMVA_
 
bool useVertices_
 
bool useVtxError_
 
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
 
edm::EDGetTokenT< reco::VertexCollectionvertices_
 
std::vector< int32_t > vtxNumber_
 vertex cuts More...
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 

Detailed Description

selects a subset of a track collection, copying extra information on demand

Author
David Lange

Definition at line 55 of file HIMultiTrackSelector.h.

Member Typedef Documentation

◆ Point

Definition at line 73 of file HIMultiTrackSelector.h.

Constructor & Destructor Documentation

◆ HIMultiTrackSelector() [1/2]

HIMultiTrackSelector::HIMultiTrackSelector ( )
explicit

constructor

Definition at line 73 of file HIMultiTrackSelector.cc.

73  {
74  useForestFromDB_ = true;
75  forest_ = nullptr;
76 }

◆ HIMultiTrackSelector() [2/2]

HIMultiTrackSelector::HIMultiTrackSelector ( const edm::ParameterSet cfg)
explicit

Definition at line 122 of file HIMultiTrackSelector.cc.

References applyAbsCutsIfNoPV_, applyAdaptedPVCuts_, applyPixelMergingCuts_, looper::cfg, chi2n_no1Dmod_par_, chi2n_par_, d0_par1_, d0_par2_, dbFileName_, dz_par1_, dz_par2_, deDxTools::esConsumes(), Exception, forest_, forestLabel_, forestToken_, forestVars_, mps_fire::i, dqmiolumiharvest::j, keepAllTracks_, max_d0_, max_d0NoPV_, max_eta_, max_lostHitFraction_, max_lostLayers_, max_minMissHitOutOrIn_, max_relpterr_, max_z0_, max_z0NoPV_, min_3Dlayers_, min_eta_, min_hits_bypass_, min_layers_, min_MVA_, min_nhits_, mvaType_, name_, nSigmaZ_, ParseForestVars(), pixel_pTMaxCut_, pixel_pTMinCut_, preFilter_, pixelTrack::qualityByName(), qualityToSet_, res_par_, setQualityBit_, AlCaHLTBitMon_QueryRunRegistry::string, useAnyMVA_, useForestFromDB_, useMVA_, useVertices_, vertexCut_, vertices_, and vtxNumber_.

123  : src_(consumes<reco::TrackCollection>(cfg.getParameter<edm::InputTag>("src"))),
124  hSrc_(consumes<TrackingRecHitCollection>(cfg.getParameter<edm::InputTag>("src"))),
125  beamspot_(consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamspot"))),
126  useVertices_(cfg.getParameter<bool>("useVertices")),
127  useVtxError_(cfg.getParameter<bool>("useVtxError"))
128 // now get the pset for each selector
129 {
130  if (useVertices_)
131  vertices_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("vertices"));
132 
133  applyPixelMergingCuts_ = false;
134  if (cfg.exists("applyPixelMergingCuts"))
135  applyPixelMergingCuts_ = cfg.getParameter<bool>("applyPixelMergingCuts");
136 
137  useAnyMVA_ = false;
138  forestLabel_ = "MVASelectorIter0";
139  std::string type = "BDTG";
140  useForestFromDB_ = true;
141  dbFileName_ = "";
142 
143  forest_ = nullptr;
144 
145  if (cfg.exists("useAnyMVA"))
146  useAnyMVA_ = cfg.getParameter<bool>("useAnyMVA");
147  if (useAnyMVA_) {
148  if (cfg.exists("mvaType"))
149  type = cfg.getParameter<std::string>("mvaType");
150  if (cfg.exists("GBRForestLabel"))
151  forestLabel_ = cfg.getParameter<std::string>("GBRForestLabel");
152  if (cfg.exists("GBRForestVars")) {
153  forestVars_ = cfg.getParameter<std::vector<std::string>>("GBRForestVars");
154  ParseForestVars();
155  }
156  if (cfg.exists("GBRForestFileName")) {
157  dbFileName_ = cfg.getParameter<std::string>("GBRForestFileName");
158  useForestFromDB_ = false;
159  }
160  mvaType_ = type;
161  }
162  if (useForestFromDB_) {
164  }
165  std::vector<edm::ParameterSet> trkSelectors(cfg.getParameter<std::vector<edm::ParameterSet>>("trackSelectors"));
166  qualityToSet_.reserve(trkSelectors.size());
167  vtxNumber_.reserve(trkSelectors.size());
168  vertexCut_.reserve(trkSelectors.size());
169  res_par_.reserve(trkSelectors.size());
170  chi2n_par_.reserve(trkSelectors.size());
171  chi2n_no1Dmod_par_.reserve(trkSelectors.size());
172  d0_par1_.reserve(trkSelectors.size());
173  dz_par1_.reserve(trkSelectors.size());
174  d0_par2_.reserve(trkSelectors.size());
175  dz_par2_.reserve(trkSelectors.size());
176  applyAdaptedPVCuts_.reserve(trkSelectors.size());
177  max_d0_.reserve(trkSelectors.size());
178  max_z0_.reserve(trkSelectors.size());
179  nSigmaZ_.reserve(trkSelectors.size());
180  pixel_pTMinCut_.reserve(trkSelectors.size());
181  pixel_pTMaxCut_.reserve(trkSelectors.size());
182  min_layers_.reserve(trkSelectors.size());
183  min_3Dlayers_.reserve(trkSelectors.size());
184  max_lostLayers_.reserve(trkSelectors.size());
185  min_hits_bypass_.reserve(trkSelectors.size());
186  applyAbsCutsIfNoPV_.reserve(trkSelectors.size());
187  max_d0NoPV_.reserve(trkSelectors.size());
188  max_z0NoPV_.reserve(trkSelectors.size());
189  preFilter_.reserve(trkSelectors.size());
190  max_relpterr_.reserve(trkSelectors.size());
191  min_nhits_.reserve(trkSelectors.size());
192  max_minMissHitOutOrIn_.reserve(trkSelectors.size());
193  max_lostHitFraction_.reserve(trkSelectors.size());
194  min_eta_.reserve(trkSelectors.size());
195  max_eta_.reserve(trkSelectors.size());
196  useMVA_.reserve(trkSelectors.size());
197  //mvaReaders_.reserve(trkSelectors.size());
198  min_MVA_.reserve(trkSelectors.size());
199  //mvaType_.reserve(trkSelectors.size());
200 
201  produces<edm::ValueMap<float>>("MVAVals");
202 
203  for (unsigned int i = 0; i < trkSelectors.size(); i++) {
204  qualityToSet_.push_back(TrackBase::undefQuality);
205  // parameters for vertex selection
206  vtxNumber_.push_back(useVertices_ ? trkSelectors[i].getParameter<int32_t>("vtxNumber") : 0);
207  vertexCut_.push_back(useVertices_ ? trkSelectors[i].getParameter<std::string>("vertexCut") : nullptr);
208  // parameters for adapted optimal cuts on chi2 and primary vertex compatibility
209  res_par_.push_back(trkSelectors[i].getParameter<std::vector<double>>("res_par"));
210  chi2n_par_.push_back(trkSelectors[i].getParameter<double>("chi2n_par"));
211  chi2n_no1Dmod_par_.push_back(trkSelectors[i].getParameter<double>("chi2n_no1Dmod_par"));
212  d0_par1_.push_back(trkSelectors[i].getParameter<std::vector<double>>("d0_par1"));
213  dz_par1_.push_back(trkSelectors[i].getParameter<std::vector<double>>("dz_par1"));
214  d0_par2_.push_back(trkSelectors[i].getParameter<std::vector<double>>("d0_par2"));
215  dz_par2_.push_back(trkSelectors[i].getParameter<std::vector<double>>("dz_par2"));
216  // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
217  applyAdaptedPVCuts_.push_back(trkSelectors[i].getParameter<bool>("applyAdaptedPVCuts"));
218  // Impact parameter absolute cuts.
219  max_d0_.push_back(trkSelectors[i].getParameter<double>("max_d0"));
220  max_z0_.push_back(trkSelectors[i].getParameter<double>("max_z0"));
221  nSigmaZ_.push_back(trkSelectors[i].getParameter<double>("nSigmaZ"));
222  // Cuts on numbers of layers with hits/3D hits/lost hits.
223  min_layers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumberLayers"));
224  min_3Dlayers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumber3DLayers"));
225  max_lostLayers_.push_back(trkSelectors[i].getParameter<uint32_t>("maxNumberLostLayers"));
226  min_hits_bypass_.push_back(trkSelectors[i].getParameter<uint32_t>("minHitsToBypassChecks"));
227  // Flag to apply absolute cuts if no PV passes the selection
228  applyAbsCutsIfNoPV_.push_back(trkSelectors[i].getParameter<bool>("applyAbsCutsIfNoPV"));
229  keepAllTracks_.push_back(trkSelectors[i].getParameter<bool>("keepAllTracks"));
230  max_relpterr_.push_back(trkSelectors[i].getParameter<double>("max_relpterr"));
231  min_nhits_.push_back(trkSelectors[i].getParameter<uint32_t>("min_nhits"));
232  max_minMissHitOutOrIn_.push_back(trkSelectors[i].existsAs<int32_t>("max_minMissHitOutOrIn")
233  ? trkSelectors[i].getParameter<int32_t>("max_minMissHitOutOrIn")
234  : 99);
235  max_lostHitFraction_.push_back(trkSelectors[i].existsAs<double>("max_lostHitFraction")
236  ? trkSelectors[i].getParameter<double>("max_lostHitFraction")
237  : 1.0);
238  min_eta_.push_back(trkSelectors[i].existsAs<double>("min_eta") ? trkSelectors[i].getParameter<double>("min_eta")
239  : -9999);
240  max_eta_.push_back(trkSelectors[i].existsAs<double>("max_eta") ? trkSelectors[i].getParameter<double>("max_eta")
241  : 9999);
242 
243  setQualityBit_.push_back(false);
244  std::string qualityStr = trkSelectors[i].getParameter<std::string>("qualityBit");
245  if (!qualityStr.empty()) {
246  setQualityBit_[i] = true;
247  qualityToSet_[i] = TrackBase::qualityByName(trkSelectors[i].getParameter<std::string>("qualityBit"));
248  }
249 
250  if (setQualityBit_[i] && (qualityToSet_[i] == TrackBase::undefQuality))
251  throw cms::Exception("Configuration")
252  << "You can't set the quality bit " << trkSelectors[i].getParameter<std::string>("qualityBit")
253  << " as it is 'undefQuality' or unknown.\n";
254 
255  if (applyAbsCutsIfNoPV_[i]) {
256  max_d0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_d0NoPV"));
257  max_z0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_z0NoPV"));
258  } else { //dummy values
259  max_d0NoPV_.push_back(0.);
260  max_z0NoPV_.push_back(0.);
261  }
262 
263  name_.push_back(trkSelectors[i].getParameter<std::string>("name"));
264 
265  preFilter_[i] = trkSelectors.size(); // no prefilter
266 
267  std::string pfName = trkSelectors[i].getParameter<std::string>("preFilterName");
268  if (!pfName.empty()) {
269  bool foundPF = false;
270  for (unsigned int j = 0; j < i; j++)
271  if (name_[j] == pfName) {
272  foundPF = true;
273  preFilter_[i] = j;
274  }
275  if (!foundPF)
276  throw cms::Exception("Configuration") << "Invalid prefilter name in HIMultiTrackSelector "
277  << trkSelectors[i].getParameter<std::string>("preFilterName");
278  }
279 
281  pixel_pTMinCut_.push_back(trkSelectors[i].getParameter<std::vector<double>>("pixel_pTMinCut"));
282  pixel_pTMaxCut_.push_back(trkSelectors[i].getParameter<std::vector<double>>("pixel_pTMaxCut"));
283  }
284 
285  // produces<std::vector<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
286  produces<edm::ValueMap<int>>(name_[i]).setBranchAlias(name_[i] + "TrackQuals");
287  if (useAnyMVA_) {
288  bool thisMVA = false;
289  if (trkSelectors[i].exists("useMVA"))
290  thisMVA = trkSelectors[i].getParameter<bool>("useMVA");
291  useMVA_.push_back(thisMVA);
292  if (thisMVA) {
293  double minVal = -1;
294  if (trkSelectors[i].exists("minMVA"))
295  minVal = trkSelectors[i].getParameter<double>("minMVA");
296  min_MVA_.push_back(minVal);
297 
298  } else {
299  min_MVA_.push_back(-9999.0);
300  }
301  } else {
302  min_MVA_.push_back(-9999.0);
303  }
304  }
305 }
std::vector< std::string > name_
std::vector< std::vector< double > > dz_par2_
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
std::vector< uint32_t > max_lostLayers_
std::vector< std::vector< double > > dz_par1_
std::vector< double > chi2n_no1Dmod_par_
std::vector< uint32_t > min_hits_bypass_
Quality qualityByName(std::string const &name)
std::vector< std::vector< double > > pixel_pTMaxCut_
std::vector< uint32_t > min_nhits_
std::vector< bool > useMVA_
std::vector< int32_t > max_minMissHitOutOrIn_
std::vector< uint32_t > min_3Dlayers_
std::vector< double > min_eta_
std::vector< double > max_eta_
std::vector< double > max_d0NoPV_
edm::EDGetTokenT< reco::VertexCollection > vertices_
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
std::vector< bool > setQualityBit_
do I have to set a quality bit?
std::vector< std::vector< double > > d0_par1_
std::vector< int32_t > vtxNumber_
vertex cuts
std::vector< double > max_z0_
edm::EDGetTokenT< TrackingRecHitCollection > hSrc_
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
std::vector< double > nSigmaZ_
std::vector< double > max_d0_
Impact parameter absolute cuts.
std::vector< bool > applyAdaptedPVCuts_
std::vector< double > max_relpterr_
std::vector< std::vector< double > > d0_par2_
std::vector< std::string > forestVars_
std::vector< int32_t > max_lostHitFraction_
std::vector< double > chi2n_par_
std::vector< std::vector< double > > res_par_
edm::EDGetTokenT< reco::BeamSpot > beamspot_
std::vector< bool > applyAbsCutsIfNoPV_
std::vector< uint32_t > min_layers_
Cuts on numbers of layers with hits/3D hits/lost hits.
std::vector< std::vector< double > > pixel_pTMinCut_
std::vector< double > max_z0NoPV_
std::vector< bool > keepAllTracks_
edm::ESGetToken< GBRForest, GBRWrapperRcd > forestToken_
std::vector< double > min_MVA_
std::vector< unsigned int > preFilter_
std::vector< reco::TrackBase::TrackQuality > qualityToSet_

◆ ~HIMultiTrackSelector()

HIMultiTrackSelector::~HIMultiTrackSelector ( )
override

destructor

Definition at line 307 of file HIMultiTrackSelector.cc.

References forest_.

307 { delete forest_; }

Member Function Documentation

◆ beginStream()

void HIMultiTrackSelector::beginStream ( edm::StreamID  )
finalprotected

Definition at line 309 of file HIMultiTrackSelector.cc.

References dbFileName_, forest_, forestLabel_, and useForestFromDB_.

309  {
310  if (!useForestFromDB_) {
311  TFile gbrfile(dbFileName_.c_str());
312  forest_ = (GBRForest *)gbrfile.Get(forestLabel_.c_str());
313  }
314 }

◆ ParseForestVars()

void HIMultiTrackSelector::ParseForestVars ( )
protected

Definition at line 78 of file HIMultiTrackSelector.cc.

References chi2n, chi2n_no1dmod, chi2perdofperlayer, dxyperdxyerror, dzperdzerror, PVValHelper::eta, etaerror, mps_fire::i, lostmidfrac, minlost, ndof, nhits, nlayers, nlayers3d, nlayerslost, relpterr, AlCaHLTBitMon_QueryRunRegistry::string, and findQualityFiles::v.

Referenced by HIMultiTrackSelector().

78  {
79  mvavars_indices.clear();
80  for (unsigned i = 0; i < forestVars_.size(); i++) {
82  int ind = -1;
83  if (v == "chi2perdofperlayer")
84  ind = chi2perdofperlayer;
85  if (v == "dxyperdxyerror")
86  ind = dxyperdxyerror;
87  if (v == "dzperdzerror")
88  ind = dzperdzerror;
89  if (v == "relpterr")
90  ind = relpterr;
91  if (v == "lostmidfrac")
92  ind = lostmidfrac;
93  if (v == "minlost")
94  ind = minlost;
95  if (v == "nhits")
96  ind = nhits;
97  if (v == "eta")
98  ind = eta;
99  if (v == "chi2n_no1dmod")
100  ind = chi2n_no1dmod;
101  if (v == "chi2n")
102  ind = chi2n;
103  if (v == "nlayerslost")
104  ind = nlayerslost;
105  if (v == "nlayers3d")
106  ind = nlayers3d;
107  if (v == "nlayers")
108  ind = nlayers;
109  if (v == "ndof")
110  ind = ndof;
111  if (v == "etaerror")
112  ind = etaerror;
113 
114  if (ind == -1)
115  edm::LogWarning("HIMultiTrackSelector")
116  << "Unknown forest variable " << v << ". Please make sure it's in the list of supported variables\n";
117 
118  mvavars_indices.push_back(ind);
119  }
120 }
std::vector< int > mvavars_indices
std::vector< std::string > forestVars_
Log< level::Warning, false > LogWarning

◆ processMVA()

void HIMultiTrackSelector::processMVA ( edm::Event evt,
const edm::EventSetup es,
std::vector< float > &  mvaVals_,
const reco::VertexCollection hVtx 
) const
protected

Definition at line 637 of file HIMultiTrackSelector.cc.

References cms::cuda::assert(), hltPixelTracks_cff::chi2, reco::TrackBase::chi2(), chi2n, chi2n_no1dmod, chi2perdofperlayer, PVValHelper::dxy, reco::TrackBase::dxy(), reco::TrackBase::dxyError(), dxyperdxyerror, PVValHelper::dz, reco::TrackBase::dz(), reco::TrackBase::dzError(), dzperdzerror, PVValHelper::eta, reco::TrackBase::eta(), etaerror, reco::TrackBase::etaError(), reco::Track::extra(), f, spr::find(), dqmMemoryStats::float, forest_, forestToken_, edm::Event::getByToken(), GBRForest::GetClassifier(), edm::EventSetup::getData(), reco::TrackBase::hitPattern(), hSrc_, mps_fire::i, lostmidfrac, SiStripPI::max, SiStripPI::min, minlost, reco::HitPattern::MISSING_INNER_HITS, reco::HitPattern::MISSING_OUTER_HITS, eostools::move(), mvavars_indices, ndof, reco::TrackBase::ndof(), nhits, nlayers, nlayers3d, nlayerslost, reco::TrackBase::normalizedChi2(), reco::TrackBase::numberOfLostHits(), reco::HitPattern::numberOfLostTrackerHits(), reco::TrackBase::numberOfValidHits(), reco::HitPattern::numberOfValidStripLayersWithMonoAndStereo(), reco::HitPattern::pixelLayersWithMeasurement(), position, reco::TrackBase::pt(), reco::TrackBase::ptError(), edm::Event::put(), reco::Track::recHitsSize(), relpterr, mathSSE::sqrt(), src_, gmt_cfi::srcTracks, reco::HitPattern::TRACK_HITS, reco::HitPattern::trackerLayersWithMeasurement(), reco::HitPattern::trackerLayersWithoutMeasurement(), useAnyMVA_, useForestFromDB_, and AlignmentTracksFromVertexSelector_cfi::vertices.

Referenced by run().

640  {
641  using namespace std;
642  using namespace edm;
643  using namespace reco;
644 
645  // Get tracks
646  Handle<TrackCollection> hSrcTrack;
647  evt.getByToken(src_, hSrcTrack);
648  const TrackCollection &srcTracks(*hSrcTrack);
649  assert(mvaVals_.size() == srcTracks.size());
650 
651  // get hits in track..
653  evt.getByToken(hSrc_, hSrcHits);
654  const TrackingRecHitCollection &srcHits(*hSrcHits);
655 
656  auto mvaValValueMap = std::make_unique<edm::ValueMap<float>>();
657  edm::ValueMap<float>::Filler mvaFiller(*mvaValValueMap);
658 
659  if (!useAnyMVA_) {
660  // mvaVals_ already initalized...
661  mvaFiller.insert(hSrcTrack, mvaVals_.begin(), mvaVals_.end());
662  mvaFiller.fill();
663  evt.put(std::move(mvaValValueMap), "MVAVals");
664  return;
665  }
666 
667  bool checkvertex =
670 
671  size_t current = 0;
672  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
673  const Track &trk = *it;
674 
675  float mvavalues[15];
676  mvavalues[ndof] = trk.ndof();
677  mvavalues[nlayers] = trk.hitPattern().trackerLayersWithMeasurement();
678  mvavalues[nlayers3d] =
681  mvavalues[chi2n_no1dmod] = trk.normalizedChi2();
682  mvavalues[chi2perdofperlayer] = mvavalues[chi2n_no1dmod] / mvavalues[nlayers];
683 
684  float chi2n1d = trk.normalizedChi2();
685  int count1dhits = 0;
686  auto ith = trk.extra()->firstRecHit();
687  auto edh = ith + trk.recHitsSize();
688  for (; ith < edh; ++ith) {
689  const TrackingRecHit &hit = srcHits[ith];
690  if (hit.dimension() == 1)
691  ++count1dhits;
692  }
693  if (count1dhits > 0) {
694  float chi2 = trk.chi2();
695  float ndof = trk.ndof();
696  chi2n1d = (chi2 + count1dhits) / float(ndof + count1dhits);
697  }
698 
699  mvavalues[chi2n] = chi2n1d; //chi2 and 1d modes
700 
701  mvavalues[eta] = trk.eta();
702  mvavalues[relpterr] = float(trk.ptError()) / std::max(float(trk.pt()), 0.000001f);
703  mvavalues[nhits] = trk.numberOfValidHits();
704 
707  mvavalues[minlost] = std::min(lostIn, lostOut);
708  mvavalues[lostmidfrac] = trk.numberOfLostHits() / (trk.numberOfValidHits() + trk.numberOfLostHits());
709 
710  mvavalues[etaerror] = trk.etaError();
711 
712  float reldz = 0;
713  float reldxy = 0;
714  if (checkvertex) {
715  int vtxind = 0; // only first vertex is taken into account for the speed purposes
716  float dxy = trk.dxy(vertices[vtxind].position()),
717  dxyE = sqrt(trk.dxyError() * trk.dxyError() + vertices[vtxind].xError() * vertices[vtxind].yError());
718  float dz = trk.dz(vertices[vtxind].position()),
719  dzE = sqrt(trk.dzError() * trk.dzError() + vertices[vtxind].zError() * vertices[vtxind].zError());
720  reldz = dz / dzE;
721  reldxy = dxy / dxyE;
722  }
723  mvavalues[dxyperdxyerror] = reldxy;
724  mvavalues[dzperdzerror] = reldz;
725 
726  std::vector<float> gbrValues;
727 
728  //fill in the gbrValues vector with the necessary variables
729  for (unsigned i = 0; i < mvavars_indices.size(); i++) {
730  gbrValues.push_back(mvavalues[mvavars_indices[i]]);
731  }
732 
733  GBRForest const *forest = forest_;
734  if (useForestFromDB_) {
735  forest = &es.getData(forestToken_);
736  }
737 
738  auto gbrVal = forest->GetClassifier(&gbrValues[0]);
739  mvaVals_[current] = gbrVal;
740  }
741  mvaFiller.insert(hSrcTrack, mvaVals_.begin(), mvaVals_.end());
742  mvaFiller.fill();
743  evt.put(std::move(mvaValValueMap), "MVAVals");
744 }
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:893
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
srcTracks
Definition: gmt_cfi.py:47
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:534
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
assert(be >=bs)
double pt() const
track transverse momentum
Definition: TrackBase.h:637
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:590
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:369
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
double dxyError() const
error on dxy
Definition: TrackBase.h:769
double dzError() const
error on dz
Definition: TrackBase.h:778
T sqrt(T t)
Definition: SSEVec.h:19
double f[11][100]
edm::EDGetTokenT< TrackingRecHitCollection > hSrc_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:801
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
std::vector< int > mvavars_indices
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
double GetClassifier(const float *vector) const
Definition: GBRForest.h:33
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:553
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
fixed size matrix
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:289
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:513
edm::ESGetToken< GBRForest, GBRWrapperRcd > forestToken_
double etaError() const
error on eta
Definition: TrackBase.h:763
def move(src, dest)
Definition: eostools.py:511
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:139

◆ produce()

void HIMultiTrackSelector::produce ( edm::Event evt,
const edm::EventSetup es 
)
inlinefinalprotected

process one event

Definition at line 75 of file HIMultiTrackSelector.h.

References writedatasetfile::run.

75 { run(evt, es); }
virtual void run(edm::Event &evt, const edm::EventSetup &es) const

◆ run()

void HIMultiTrackSelector::run ( edm::Event evt,
const edm::EventSetup es 
) const
protectedvirtual

Definition at line 316 of file HIMultiTrackSelector.cc.

References beamspot_, trigObjTnPSource_cfi::filler, edm::Event::getByToken(), muons_cff::highPurity, hSrc_, mps_fire::i, dqmiolumiharvest::j, keepAllTracks_, LogTrace, loose, eostools::move(), name_, convertSQLiteXML::ok, hiPixelPairStep_cff::points, preFilter_, processMVA(), reco::TrackBase::pt(), edm::Event::put(), reco::TrackBase::qualityMask(), qualityToSet_, select(), selectVertices(), setQualityBit_, src_, gmt_cfi::srcTracks, useAnyMVA_, useVertices_, and vertices_.

316  {
317  using namespace std;
318  using namespace edm;
319  using namespace reco;
320 
321  // Get tracks
322  Handle<TrackCollection> hSrcTrack;
323  evt.getByToken(src_, hSrcTrack);
324  const TrackCollection &srcTracks(*hSrcTrack);
325 
326  // get hits in track..
328  evt.getByToken(hSrc_, hSrcHits);
329  const TrackingRecHitCollection &srcHits(*hSrcHits);
330 
331  // looking for the beam spot
333  evt.getByToken(beamspot_, hBsp);
334  const reco::BeamSpot &vertexBeamSpot(*hBsp);
335 
336  // Select good primary vertices for use in subsequent track selection
338  if (useVertices_)
339  evt.getByToken(vertices_, hVtx);
340 
341  unsigned int trkSize = srcTracks.size();
342  std::vector<int> selTracksSave(qualityToSet_.size() * trkSize, 0);
343 
344  std::vector<float> mvaVals_(srcTracks.size(), -99.f);
345  processMVA(evt, es, mvaVals_, *hVtx);
346 
347  for (unsigned int i = 0; i < qualityToSet_.size(); i++) {
348  std::vector<int> selTracks(trkSize, 0);
349  auto selTracksValueMap = std::make_unique<edm::ValueMap<int>>();
350  edm::ValueMap<int>::Filler filler(*selTracksValueMap);
351 
352  std::vector<Point> points;
353  std::vector<float> vterr, vzerr;
354  if (useVertices_)
355  selectVertices(i, *hVtx, points, vterr, vzerr);
356 
357  // Loop over tracks
358  size_t current = 0;
359  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
360  const Track &trk = *it;
361  // Check if this track passes cuts
362 
363  LogTrace("TrackSelection") << "ready to check track with pt=" << trk.pt();
364 
365  //already removed
366  bool ok = true;
367  float mvaVal = 0;
368  if (preFilter_[i] < i && selTracksSave[preFilter_[i] * trkSize + current] < 0) {
369  selTracks[current] = -1;
370  ok = false;
371  if (!keepAllTracks_[i])
372  continue;
373  } else {
374  if (useAnyMVA_)
375  mvaVal = mvaVals_[current];
376  ok = select(i, vertexBeamSpot, srcHits, trk, points, vterr, vzerr, mvaVal);
377  if (!ok) {
378  LogTrace("TrackSelection") << "track with pt=" << trk.pt() << " NOT selected";
379  if (!keepAllTracks_[i]) {
380  selTracks[current] = -1;
381  continue;
382  }
383  } else
384  LogTrace("TrackSelection") << "track with pt=" << trk.pt() << " selected";
385  }
386 
387  if (preFilter_[i] < i) {
388  selTracks[current] = selTracksSave[preFilter_[i] * trkSize + current];
389  } else {
390  selTracks[current] = trk.qualityMask();
391  }
392  if (ok && setQualityBit_[i]) {
393  selTracks[current] = (selTracks[current] | (1 << qualityToSet_[i]));
394  if (qualityToSet_[i] == TrackBase::tight) {
395  selTracks[current] = (selTracks[current] | (1 << TrackBase::loose));
396  } else if (qualityToSet_[i] == TrackBase::highPurity) {
397  selTracks[current] = (selTracks[current] | (1 << TrackBase::loose));
398  selTracks[current] = (selTracks[current] | (1 << TrackBase::tight));
399  }
400 
401  if (!points.empty()) {
402  if (qualityToSet_[i] == TrackBase::loose) {
403  selTracks[current] = (selTracks[current] | (1 << TrackBase::looseSetWithPV));
404  } else if (qualityToSet_[i] == TrackBase::highPurity) {
405  selTracks[current] = (selTracks[current] | (1 << TrackBase::looseSetWithPV));
406  selTracks[current] = (selTracks[current] | (1 << TrackBase::highPuritySetWithPV));
407  }
408  }
409  }
410  }
411  for (unsigned int j = 0; j < trkSize; j++)
412  selTracksSave[j + i * trkSize] = selTracks[j];
413  filler.insert(hSrcTrack, selTracks.begin(), selTracks.end());
414  filler.fill();
415 
416  // evt.put(std::move(selTracks),name_[i]);
417  evt.put(std::move(selTracksValueMap), name_[i]);
418  }
419 }
std::vector< std::string > name_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
srcTracks
Definition: gmt_cfi.py:47
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define LogTrace(id)
double pt() const
track transverse momentum
Definition: TrackBase.h:637
void selectVertices(unsigned int tsNum, const reco::VertexCollection &vtxs, std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr) const
bool select(unsigned tsNum, const reco::BeamSpot &vertexBeamSpot, const TrackingRecHitCollection &recHits, const reco::Track &tk, const std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr, double mvaVal) const
return class, or -1 if rejected
edm::EDGetTokenT< reco::VertexCollection > vertices_
std::vector< bool > setQualityBit_
do I have to set a quality bit?
edm::EDGetTokenT< TrackingRecHitCollection > hSrc_
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
int qualityMask() const
Definition: TrackBase.h:843
edm::EDGetTokenT< reco::BeamSpot > beamspot_
fixed size matrix
HLT enums.
constexpr auto loose
std::vector< bool > keepAllTracks_
std::vector< unsigned int > preFilter_
std::vector< reco::TrackBase::TrackQuality > qualityToSet_
void processMVA(edm::Event &evt, const edm::EventSetup &es, std::vector< float > &mvaVals_, const reco::VertexCollection &hVtx) const
def move(src, dest)
Definition: eostools.py:511

◆ select()

bool HIMultiTrackSelector::select ( unsigned  tsNum,
const reco::BeamSpot vertexBeamSpot,
const TrackingRecHitCollection recHits,
const reco::Track tk,
const std::vector< Point > &  points,
std::vector< float > &  vterr,
std::vector< float > &  vzerr,
double  mvaVal 
) const
protected

return class, or -1 if rejected

Definition at line 421 of file HIMultiTrackSelector.cc.

References funct::abs(), applyAbsCutsIfNoPV_, applyAdaptedPVCuts_, applyPixelMergingCuts_, hltPixelTracks_cff::chi2, reco::TrackBase::chi2(), chi2n, chi2n_no1Dmod_par_, chi2n_par_, d0, d0_par1_, d0_par2_, TrackSplittingMonitor_cfi::d0Cut, reco::TrackBase::d0Error(), reco::TrackBase::dxy(), PVValHelper::dz, reco::TrackBase::dz(), dz_par1_, dz_par2_, TrackSplittingMonitor_cfi::dzCut, reco::TrackBase::dzError(), mps_fire::end, PVValHelper::eta, reco::TrackBase::eta(), reco::Track::extra(), f, dqmMemoryStats::float, reco::TrackBase::hitPattern(), gpuVertexFinder::iv, LogDebug, LogTrace, SiStripPI::max, max_d0_, max_d0NoPV_, max_eta_, max_lostHitFraction_, max_lostLayers_, max_minMissHitOutOrIn_, max_relpterr_, max_z0_, max_z0NoPV_, SiStripPI::min, min_3Dlayers_, min_eta_, min_hits_bypass_, min_layers_, min_MVA_, min_nhits_, reco::HitPattern::MISSING_INNER_HITS, reco::HitPattern::MISSING_OUTER_HITS, ndof, reco::TrackBase::ndof(), nhits, nlayers, reco::TrackBase::normalizedChi2(), nSigmaZ_, reco::TrackBase::numberOfLostHits(), reco::HitPattern::numberOfLostTrackerHits(), reco::TrackBase::numberOfValidHits(), reco::HitPattern::numberOfValidStripLayersWithMonoAndStereo(), pixel_pTMaxCut_, pixel_pTMinCut_, reco::HitPattern::pixelLayersWithMeasurement(), point, hiPixelPairStep_cff::points, reco::BeamSpot::position(), funct::pow(), DiDispStaMuonMonitor_cfi::pt, reco::TrackBase::pt(), reco::TrackBase::ptError(), FastTrackerRecHitMaskProducer_cfi::recHits, reco::Track::recHitsSize(), relpterr, res_par_, reco::BeamSpot::sigmaZ(), mathSSE::sqrt(), reco::HitPattern::TRACK_HITS, reco::HitPattern::trackerLayersWithMeasurement(), reco::HitPattern::trackerLayersWithoutMeasurement(), useAnyMVA_, useMVA_, and useVtxError_.

Referenced by run().

428  {
429  // Decide if the given track passes selection cuts.
430 
431  using namespace std;
432 
433  //cuts on number of valid hits
434  auto nhits = tk.numberOfValidHits();
435  if (nhits >= min_hits_bypass_[tsNum])
436  return true;
437  if (nhits < min_nhits_[tsNum])
438  return false;
439 
440  if (tk.ndof() < 1E-5)
441  return false;
442 
444  //Adding the MVA selection before any other cut//
446  if (useAnyMVA_ && useMVA_[tsNum]) {
447  if (mvaVal < min_MVA_[tsNum])
448  return false;
449  else
450  return true;
451  }
453  //End of MVA selection section//
455 
456  // Cuts on numbers of layers with hits/3D hits/lost hits.
458  uint32_t nlayers3D =
461  LogDebug("TrackSelection") << "cuts on nlayers: " << nlayers << " " << nlayers3D << " " << nlayersLost << " vs "
462  << min_layers_[tsNum] << " " << min_3Dlayers_[tsNum] << " " << max_lostLayers_[tsNum];
463  if (nlayers < min_layers_[tsNum])
464  return false;
465  if (nlayers3D < min_3Dlayers_[tsNum])
466  return false;
467  if (nlayersLost > max_lostLayers_[tsNum])
468  return false;
469  LogTrace("TrackSelection") << "cuts on nlayers passed";
470 
471  float chi2n = tk.normalizedChi2();
472  float chi2n_no1Dmod = chi2n;
473 
474  int count1dhits = 0;
475  auto ith = tk.extra()->firstRecHit();
476  auto edh = ith + tk.recHitsSize();
477  for (; ith < edh; ++ith) {
478  const TrackingRecHit &hit = recHits[ith];
479  if (hit.dimension() == 1)
480  ++count1dhits;
481  }
482  if (count1dhits > 0) {
483  float chi2 = tk.chi2();
484  float ndof = tk.ndof();
485  chi2n = (chi2 + count1dhits) / float(ndof + count1dhits);
486  }
487  // For each 1D rechit, the chi^2 and ndof is increased by one. This is a way of retaining approximately
488  // the same normalized chi^2 distribution as with 2D rechits.
489  if (chi2n > chi2n_par_[tsNum] * nlayers)
490  return false;
491 
492  if (chi2n_no1Dmod > chi2n_no1Dmod_par_[tsNum] * nlayers)
493  return false;
494 
495  // Get track parameters
496  float pt = std::max(float(tk.pt()), 0.000001f);
497  float eta = tk.eta();
498  if (eta < min_eta_[tsNum] || eta > max_eta_[tsNum])
499  return false;
500 
501  //cuts on relative error on pt
502  float relpterr = float(tk.ptError()) / pt;
503  if (relpterr > max_relpterr_[tsNum])
504  return false;
505 
508  int minLost = std::min(lostIn, lostOut);
509  if (minLost > max_minMissHitOutOrIn_[tsNum])
510  return false;
511  float lostMidFrac =
512  tk.numberOfLostHits() == 0 ? 0. : tk.numberOfLostHits() / (tk.numberOfValidHits() + tk.numberOfLostHits());
513  if (lostMidFrac > max_lostHitFraction_[tsNum])
514  return false;
515 
516  // Pixel Track Merging pT dependent cuts
518  // hard cut at absolute min/max pt
519  if (pt < pixel_pTMinCut_[tsNum][0])
520  return false;
521  if (pt > pixel_pTMaxCut_[tsNum][0])
522  return false;
523  // tapering cuts with chi2n_no1Dmod
524  double pTMaxCutPos = (pixel_pTMaxCut_[tsNum][0] - pt) / (pixel_pTMaxCut_[tsNum][0] - pixel_pTMaxCut_[tsNum][1]);
525  double pTMinCutPos = (pt - pixel_pTMinCut_[tsNum][0]) / (pixel_pTMinCut_[tsNum][1] - pixel_pTMinCut_[tsNum][0]);
526  if (pt > pixel_pTMaxCut_[tsNum][1] &&
527  chi2n_no1Dmod > pixel_pTMaxCut_[tsNum][2] * nlayers * pow(pTMaxCutPos, pixel_pTMaxCut_[tsNum][3]))
528  return false;
529  if (pt < pixel_pTMinCut_[tsNum][1] &&
530  chi2n_no1Dmod > pixel_pTMinCut_[tsNum][2] * nlayers * pow(pTMinCutPos, pixel_pTMinCut_[tsNum][3]))
531  return false;
532  }
533 
534  //other track parameters
535  float d0 = -tk.dxy(vertexBeamSpot.position()), d0E = tk.d0Error(), dz = tk.dz(vertexBeamSpot.position()),
536  dzE = tk.dzError();
537 
538  // parametrized d0 resolution for the track pt
539  float nomd0E = sqrt(res_par_[tsNum][0] * res_par_[tsNum][0] + (res_par_[tsNum][1] / pt) * (res_par_[tsNum][1] / pt));
540  // parametrized z0 resolution for the track pt and eta
541  float nomdzE = nomd0E * (std::cosh(eta));
542 
543  float dzCut = std::min(powN(dz_par1_[tsNum][0] * nlayers, int(dz_par1_[tsNum][1] + 0.5)) * nomdzE,
544  powN(dz_par2_[tsNum][0] * nlayers, int(dz_par2_[tsNum][1] + 0.5)) * dzE);
545  float d0Cut = std::min(powN(d0_par1_[tsNum][0] * nlayers, int(d0_par1_[tsNum][1] + 0.5)) * nomd0E,
546  powN(d0_par2_[tsNum][0] * nlayers, int(d0_par2_[tsNum][1] + 0.5)) * d0E);
547 
548  // ---- PrimaryVertex compatibility cut
549  bool primaryVertexZCompatibility(false);
550  bool primaryVertexD0Compatibility(false);
551 
552  if (points.empty()) { //If not primaryVertices are reconstructed, check just the compatibility with the BS
553  //z0 within (n sigma + dzCut) of the beam spot z, if no good vertex is found
554  if (abs(dz) < hypot(vertexBeamSpot.sigmaZ() * nSigmaZ_[tsNum], dzCut))
555  primaryVertexZCompatibility = true;
556  // d0 compatibility with beam line
557  if (abs(d0) < d0Cut)
558  primaryVertexD0Compatibility = true;
559  }
560 
561  int iv = 0;
562  for (std::vector<Point>::const_iterator point = points.begin(), end = points.end(); point != end; ++point) {
563  LogTrace("TrackSelection") << "Test track w.r.t. vertex with z position " << point->z();
564  if (primaryVertexZCompatibility && primaryVertexD0Compatibility)
565  break;
566  float dzPV = tk.dz(*point); //re-evaluate the dz with respect to the vertex position
567  float d0PV = tk.dxy(*point); //re-evaluate the dxy with respect to the vertex position
568  if (useVtxError_) {
569  float dzErrPV = std::sqrt(dzE * dzE + vzerr[iv] * vzerr[iv]); // include vertex error in z
570  float d0ErrPV = std::sqrt(d0E * d0E + vterr[iv] * vterr[iv]); // include vertex error in xy
571  iv++;
572  if (abs(dzPV) < dz_par1_[tsNum][0] * pow(nlayers, dz_par1_[tsNum][1]) * nomdzE &&
573  abs(dzPV) < dz_par2_[tsNum][0] * pow(nlayers, dz_par2_[tsNum][1]) * dzErrPV && abs(dzPV) < max_z0_[tsNum])
574  primaryVertexZCompatibility = true;
575  if (abs(d0PV) < d0_par1_[tsNum][0] * pow(nlayers, d0_par1_[tsNum][1]) * nomd0E &&
576  abs(d0PV) < d0_par2_[tsNum][0] * pow(nlayers, d0_par2_[tsNum][1]) * d0ErrPV && abs(d0PV) < max_d0_[tsNum])
577  primaryVertexD0Compatibility = true;
578  } else {
579  if (abs(dzPV) < dzCut)
580  primaryVertexZCompatibility = true;
581  if (abs(d0PV) < d0Cut)
582  primaryVertexD0Compatibility = true;
583  }
584  LogTrace("TrackSelection") << "distances " << dzPV << " " << d0PV << " vs " << dzCut << " " << d0Cut;
585  }
586 
587  if (points.empty() && applyAbsCutsIfNoPV_[tsNum]) {
588  if (abs(dz) > max_z0NoPV_[tsNum] || abs(d0) > max_d0NoPV_[tsNum])
589  return false;
590  } else {
591  // Absolute cuts on all tracks impact parameters with respect to beam-spot.
592  // If BS is not compatible, verify if at least the reco-vertex is compatible (useful for incorrect BS settings)
593  if (abs(d0) > max_d0_[tsNum] && !primaryVertexD0Compatibility)
594  return false;
595  LogTrace("TrackSelection") << "absolute cuts on d0 passed";
596  if (abs(dz) > max_z0_[tsNum] && !primaryVertexZCompatibility)
597  return false;
598  LogTrace("TrackSelection") << "absolute cuts on dz passed";
599  }
600 
601  LogTrace("TrackSelection") << "cuts on PV: apply adapted PV cuts? " << applyAdaptedPVCuts_[tsNum]
602  << " d0 compatibility? " << primaryVertexD0Compatibility << " z compatibility? "
603  << primaryVertexZCompatibility;
604 
605  if (applyAdaptedPVCuts_[tsNum]) {
606  return (primaryVertexD0Compatibility && primaryVertexZCompatibility);
607  } else {
608  return true;
609  }
610 }
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:893
std::vector< std::vector< double > > dz_par2_
std::vector< uint32_t > max_lostLayers_
std::vector< std::vector< double > > dz_par1_
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
int32_t *__restrict__ iv
std::vector< double > chi2n_no1Dmod_par_
std::vector< uint32_t > min_hits_bypass_
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
const Point & position() const
position
Definition: BeamSpot.h:59
std::vector< std::vector< double > > pixel_pTMaxCut_
std::vector< uint32_t > min_nhits_
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:534
std::vector< bool > useMVA_
std::vector< int32_t > max_minMissHitOutOrIn_
std::vector< uint32_t > min_3Dlayers_
#define LogTrace(id)
std::vector< double > min_eta_
std::vector< double > max_eta_
double pt() const
track transverse momentum
Definition: TrackBase.h:637
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:590
std::vector< double > max_d0NoPV_
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:369
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
double dzError() const
error on dz
Definition: TrackBase.h:778
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< std::vector< double > > d0_par1_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > max_z0_
double f[11][100]
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:801
std::vector< double > nSigmaZ_
std::vector< double > max_d0_
Impact parameter absolute cuts.
std::vector< bool > applyAdaptedPVCuts_
std::vector< double > max_relpterr_
std::vector< std::vector< double > > d0_par2_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
static constexpr float d0
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:553
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
double sigmaZ() const
sigma z
Definition: BeamSpot.h:76
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
std::vector< int32_t > max_lostHitFraction_
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:593
std::vector< double > chi2n_par_
std::vector< std::vector< double > > res_par_
std::vector< bool > applyAbsCutsIfNoPV_
std::vector< uint32_t > min_layers_
Cuts on numbers of layers with hits/3D hits/lost hits.
std::vector< std::vector< double > > pixel_pTMinCut_
std::vector< double > max_z0NoPV_
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:513
double d0Error() const
error on d0
Definition: TrackBase.h:772
std::vector< double > min_MVA_
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
#define LogDebug(id)
const TrackExtraRef & extra() const
reference to "extra" object
Definition: Track.h:139

◆ selectVertices()

void HIMultiTrackSelector::selectVertices ( unsigned int  tsNum,
const reco::VertexCollection vtxs,
std::vector< Point > &  points,
std::vector< float > &  vterr,
std::vector< float > &  vzerr 
) const
protected

Definition at line 612 of file HIMultiTrackSelector.cc.

References LogDebug, LogTrace, hiPixelPairStep_cff::points, mathSSE::sqrt(), vertexCut_, extraflags_cff::vtx, and vtxNumber_.

Referenced by run().

616  {
617  // Select good primary vertices
618  using namespace reco;
619  int32_t toTake = vtxNumber_[tsNum];
620  for (VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
621  LogDebug("SelectVertex") << " select vertex with z position " << it->z() << " " << it->chi2() << " " << it->ndof()
622  << " " << TMath::Prob(it->chi2(), static_cast<int32_t>(it->ndof()));
623  Vertex vtx = *it;
624  bool pass = vertexCut_[tsNum](vtx);
625  if (pass) {
626  points.push_back(it->position());
627  vterr.push_back(sqrt(it->yError() * it->xError()));
628  vzerr.push_back(it->zError());
629  LogTrace("SelectVertex") << " SELECTED vertex with z position " << it->z();
630  toTake--;
631  if (toTake == 0)
632  break;
633  }
634  }
635 }
#define LogTrace(id)
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< int32_t > vtxNumber_
vertex cuts
fixed size matrix
#define LogDebug(id)

Member Data Documentation

◆ applyAbsCutsIfNoPV_

std::vector<bool> HIMultiTrackSelector::applyAbsCutsIfNoPV_
protected

Definition at line 159 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ applyAdaptedPVCuts_

std::vector<bool> HIMultiTrackSelector::applyAdaptedPVCuts_
protected

Definition at line 129 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ applyPixelMergingCuts_

bool HIMultiTrackSelector::applyPixelMergingCuts_
protected

Definition at line 109 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ beamspot_

edm::EDGetTokenT<reco::BeamSpot> HIMultiTrackSelector::beamspot_
protected

Definition at line 102 of file HIMultiTrackSelector.h.

Referenced by run().

◆ chi2n_no1Dmod_par_

std::vector<double> HIMultiTrackSelector::chi2n_no1Dmod_par_
protected

Definition at line 123 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ chi2n_par_

std::vector<double> HIMultiTrackSelector::chi2n_par_
protected

Definition at line 122 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ d0_par1_

std::vector<std::vector<double> > HIMultiTrackSelector::d0_par1_
protected

Definition at line 124 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ d0_par2_

std::vector<std::vector<double> > HIMultiTrackSelector::d0_par2_
protected

Definition at line 126 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ dbFileName_

std::string HIMultiTrackSelector::dbFileName_
protected

Definition at line 182 of file HIMultiTrackSelector.h.

Referenced by beginStream(), and HIMultiTrackSelector().

◆ dz_par1_

std::vector<std::vector<double> > HIMultiTrackSelector::dz_par1_
protected

Definition at line 125 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ dz_par2_

std::vector<std::vector<double> > HIMultiTrackSelector::dz_par2_
protected

Definition at line 127 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ forest_

GBRForest* HIMultiTrackSelector::forest_
protected

◆ forestLabel_

std::string HIMultiTrackSelector::forestLabel_
protected

Definition at line 177 of file HIMultiTrackSelector.h.

Referenced by beginStream(), and HIMultiTrackSelector().

◆ forestToken_

edm::ESGetToken<GBRForest, GBRWrapperRcd> HIMultiTrackSelector::forestToken_
protected

Definition at line 180 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and processMVA().

◆ forestVars_

std::vector<std::string> HIMultiTrackSelector::forestVars_
protected

Definition at line 178 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector().

◆ hSrc_

edm::EDGetTokenT<TrackingRecHitCollection> HIMultiTrackSelector::hSrc_
protected

Definition at line 101 of file HIMultiTrackSelector.h.

Referenced by processMVA(), and run().

◆ keepAllTracks_

std::vector<bool> HIMultiTrackSelector::keepAllTracks_
protected

Definition at line 161 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ max_d0_

std::vector<double> HIMultiTrackSelector::max_d0_
protected

Impact parameter absolute cuts.

Definition at line 132 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_d0NoPV_

std::vector<double> HIMultiTrackSelector::max_d0NoPV_
protected

Definition at line 157 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_eta_

std::vector<double> HIMultiTrackSelector::max_eta_
protected

Definition at line 154 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_lostHitFraction_

std::vector<int32_t> HIMultiTrackSelector::max_lostHitFraction_
protected

Definition at line 151 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_lostLayers_

std::vector<uint32_t> HIMultiTrackSelector::max_lostLayers_
protected

Definition at line 143 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_minMissHitOutOrIn_

std::vector<int32_t> HIMultiTrackSelector::max_minMissHitOutOrIn_
protected

Definition at line 150 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_relpterr_

std::vector<double> HIMultiTrackSelector::max_relpterr_
protected

Definition at line 147 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_z0_

std::vector<double> HIMultiTrackSelector::max_z0_
protected

Definition at line 133 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ max_z0NoPV_

std::vector<double> HIMultiTrackSelector::max_z0NoPV_
protected

Definition at line 158 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ min_3Dlayers_

std::vector<uint32_t> HIMultiTrackSelector::min_3Dlayers_
protected

Definition at line 142 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ min_eta_

std::vector<double> HIMultiTrackSelector::min_eta_
protected

Definition at line 153 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ min_hits_bypass_

std::vector<uint32_t> HIMultiTrackSelector::min_hits_bypass_
protected

Definition at line 144 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ min_layers_

std::vector<uint32_t> HIMultiTrackSelector::min_layers_
protected

Cuts on numbers of layers with hits/3D hits/lost hits.

Definition at line 141 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ min_MVA_

std::vector<double> HIMultiTrackSelector::min_MVA_
protected

Definition at line 173 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ min_nhits_

std::vector<uint32_t> HIMultiTrackSelector::min_nhits_
protected

Definition at line 148 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ mvaType_

std::string HIMultiTrackSelector::mvaType_
protected

Definition at line 176 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector().

◆ mvavars_indices

std::vector<int> HIMultiTrackSelector::mvavars_indices
protected

Definition at line 171 of file HIMultiTrackSelector.h.

Referenced by processMVA().

◆ name_

std::vector<std::string> HIMultiTrackSelector::name_
protected

Definition at line 165 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ nSigmaZ_

std::vector<double> HIMultiTrackSelector::nSigmaZ_
protected

Definition at line 134 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ pixel_pTMaxCut_

std::vector<std::vector<double> > HIMultiTrackSelector::pixel_pTMaxCut_
protected

Definition at line 138 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ pixel_pTMinCut_

std::vector<std::vector<double> > HIMultiTrackSelector::pixel_pTMinCut_
protected

Definition at line 137 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ preFilter_

std::vector<unsigned int> HIMultiTrackSelector::preFilter_
protected

Definition at line 164 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ qualityToSet_

std::vector<reco::TrackBase::TrackQuality> HIMultiTrackSelector::qualityToSet_
protected

Definition at line 113 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ res_par_

std::vector<std::vector<double> > HIMultiTrackSelector::res_par_
protected

Definition at line 121 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ setQualityBit_

std::vector<bool> HIMultiTrackSelector::setQualityBit_
protected

do I have to set a quality bit?

Definition at line 112 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ src_

edm::EDGetTokenT<reco::TrackCollection> HIMultiTrackSelector::src_
protected

source collection label

Definition at line 100 of file HIMultiTrackSelector.h.

Referenced by processMVA(), and run().

◆ useAnyMVA_

bool HIMultiTrackSelector::useAnyMVA_
protected

Definition at line 105 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), processMVA(), run(), and select().

◆ useForestFromDB_

bool HIMultiTrackSelector::useForestFromDB_
protected

Definition at line 181 of file HIMultiTrackSelector.h.

Referenced by beginStream(), HIMultiTrackSelector(), and processMVA().

◆ useMVA_

std::vector<bool> HIMultiTrackSelector::useMVA_
protected

Definition at line 168 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and select().

◆ useVertices_

bool HIMultiTrackSelector::useVertices_
protected

Definition at line 103 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ useVtxError_

bool HIMultiTrackSelector::useVtxError_
protected

Definition at line 104 of file HIMultiTrackSelector.h.

Referenced by select().

◆ vertexCut_

std::vector<StringCutObjectSelector<reco::Vertex> > HIMultiTrackSelector::vertexCut_
protected

Definition at line 118 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and selectVertices().

◆ vertices_

edm::EDGetTokenT<reco::VertexCollection> HIMultiTrackSelector::vertices_
protected

Definition at line 106 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and run().

◆ vtxNumber_

std::vector<int32_t> HIMultiTrackSelector::vtxNumber_
protected

vertex cuts

Definition at line 116 of file HIMultiTrackSelector.h.

Referenced by HIMultiTrackSelector(), and selectVertices().