CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
List of all members | Public Types | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes
MultiTrackSelector Class Reference

#include <MultiTrackSelector.h>

Inheritance diagram for MultiTrackSelector:
edm::stream::EDProducer<> AnalyticalTrackSelector

Public Types

using MVACollection = std::vector< float >
 
using QualityMaskCollection = std::vector< unsigned char >
 
- 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
 

Public Member Functions

 MultiTrackSelector ()
 constructor More...
 
 MultiTrackSelector (const edm::ParameterSet &cfg)
 
 ~MultiTrackSelector () 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
 
Point getBestVertex (const reco::TrackBaseRef, const reco::VertexCollection) const
 
void processMVA (edm::Event &evt, const edm::EventSetup &es, const reco::BeamSpot &beamspot, const reco::VertexCollection &vertices, int selIndex, std::vector< float > &mvaVals_, bool writeIt=false) 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_
 
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_
 
std::vector< GBRForest * > forest_
 
std::vector< std::string > forestLabel_
 
std::vector< edm::ESGetToken
< GBRForest, GBRWrapperRcd > > 
forestToken_
 
edm::EDGetTokenT
< TrackingRecHitCollection
hSrc_
 
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::vector< std::string > mvaType_
 
std::vector< std::string > name_
 
std::vector< double > nSigmaZ_
 
std::vector< unsigned int > preFilter_
 
std::vector
< reco::TrackBase::TrackQuality
qualityToSet_
 
std::vector< std::vector
< double > > 
res_par_
 
std::vector< bool > setQualityBit_
 do I have to set a quality bit? More...
 
edm::EDGetTokenT
< reco::TrackCollection
src_
 source collection label More...
 
bool useAnyMVA_
 
bool useForestFromDB_
 
std::vector< bool > useMVA_
 
std::vector< bool > useMVAonly_
 
bool useVertices_
 
bool useVtxError_
 
std::vector
< StringCutObjectSelector
< reco::Vertex > > 
vertexCut_
 
edm::EDGetTokenT
< reco::VertexCollection
vertices_
 
std::vector< int32_t > vtxNumber_
 vertex cuts More...
 

Detailed Description

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

Author
David Lange

Definition at line 36 of file MultiTrackSelector.h.

Member Typedef Documentation

using MultiTrackSelector::MVACollection = std::vector<float>

Definition at line 45 of file MultiTrackSelector.h.

Definition at line 57 of file MultiTrackSelector.h.

using MultiTrackSelector::QualityMaskCollection = std::vector<unsigned char>

Definition at line 46 of file MultiTrackSelector.h.

Constructor & Destructor Documentation

MultiTrackSelector::MultiTrackSelector ( )
explicit

constructor

Definition at line 16 of file MultiTrackSelector.cc.

16 { useForestFromDB_ = true; }
MultiTrackSelector::MultiTrackSelector ( const edm::ParameterSet cfg)
explicit

Definition at line 18 of file MultiTrackSelector.cc.

References applyAbsCutsIfNoPV_, applyAdaptedPVCuts_, chi2n_no1Dmod_par_, chi2n_par_, d0_par1_, d0_par2_, dbFileName_, dz_par1_, dz_par2_, Exception, edm::ParameterSet::exists(), forest_, forestLabel_, forestToken_, edm::ParameterSet::getParameter(), mps_fire::i, dqmiolumiharvest::j, keepAllTracks_, label, 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_, preFilter_, pixelTrack::qualityByName(), qualityToSet_, res_par_, setQualityBit_, AlCaHLTBitMon_QueryRunRegistry::string, useAnyMVA_, useForestFromDB_, useMVA_, useMVAonly_, useVertices_, useVtxError_, edm::vector_transform(), vertexCut_, vertices_, and vtxNumber_.

19  : src_(consumes<reco::TrackCollection>(cfg.getParameter<edm::InputTag>("src"))),
20  hSrc_(consumes<TrackingRecHitCollection>(cfg.getParameter<edm::InputTag>("src"))),
21  beamspot_(consumes<reco::BeamSpot>(cfg.getParameter<edm::InputTag>("beamspot"))),
22  useVertices_(cfg.getParameter<bool>("useVertices")),
23  useVtxError_(cfg.getParameter<bool>("useVtxError"))
24 // now get the pset for each selector
25 {
26  if (useVertices_)
27  vertices_ = consumes<reco::VertexCollection>(cfg.getParameter<edm::InputTag>("vertices"));
28  if (useVtxError_) {
29  edm::LogWarning("MultiTRackSelector") << "you are executing buggy code, if intentional please help to fix it";
30  }
31  useAnyMVA_ = false;
32  useForestFromDB_ = true;
33  dbFileName_ = "";
34 
35  if (cfg.exists("useAnyMVA"))
36  useAnyMVA_ = cfg.getParameter<bool>("useAnyMVA");
37 
38  if (useAnyMVA_) {
39  if (cfg.exists("GBRForestFileName")) {
40  dbFileName_ = cfg.getParameter<std::string>("GBRForestFileName");
41  useForestFromDB_ = false;
42  }
43  }
44  std::vector<edm::ParameterSet> trkSelectors(cfg.getParameter<std::vector<edm::ParameterSet>>("trackSelectors"));
45  qualityToSet_.reserve(trkSelectors.size());
46  vtxNumber_.reserve(trkSelectors.size());
47  vertexCut_.reserve(trkSelectors.size());
48  res_par_.reserve(trkSelectors.size());
49  chi2n_par_.reserve(trkSelectors.size());
50  chi2n_no1Dmod_par_.reserve(trkSelectors.size());
51  d0_par1_.reserve(trkSelectors.size());
52  dz_par1_.reserve(trkSelectors.size());
53  d0_par2_.reserve(trkSelectors.size());
54  dz_par2_.reserve(trkSelectors.size());
55  applyAdaptedPVCuts_.reserve(trkSelectors.size());
56  max_d0_.reserve(trkSelectors.size());
57  max_z0_.reserve(trkSelectors.size());
58  nSigmaZ_.reserve(trkSelectors.size());
59  min_layers_.reserve(trkSelectors.size());
60  min_3Dlayers_.reserve(trkSelectors.size());
61  max_lostLayers_.reserve(trkSelectors.size());
62  min_hits_bypass_.reserve(trkSelectors.size());
63  applyAbsCutsIfNoPV_.reserve(trkSelectors.size());
64  max_d0NoPV_.reserve(trkSelectors.size());
65  max_z0NoPV_.reserve(trkSelectors.size());
66  preFilter_.reserve(trkSelectors.size());
67  max_relpterr_.reserve(trkSelectors.size());
68  min_nhits_.reserve(trkSelectors.size());
69  max_minMissHitOutOrIn_.reserve(trkSelectors.size());
70  max_lostHitFraction_.reserve(trkSelectors.size());
71  min_eta_.reserve(trkSelectors.size());
72  max_eta_.reserve(trkSelectors.size());
73  useMVA_.reserve(trkSelectors.size());
74  useMVAonly_.reserve(trkSelectors.size());
75  //mvaReaders_.reserve(trkSelectors.size());
76  min_MVA_.reserve(trkSelectors.size());
77  mvaType_.reserve(trkSelectors.size());
78  forestLabel_.reserve(trkSelectors.size());
79  forest_.reserve(trkSelectors.size());
80 
81  produces<edm::ValueMap<float>>("MVAVals");
82 
83  //foward compatibility
84  produces<MVACollection>("MVAValues");
85 
86  for (unsigned int i = 0; i < trkSelectors.size(); i++) {
87  qualityToSet_.push_back(TrackBase::undefQuality);
88  // parameters for vertex selection
89  vtxNumber_.push_back(useVertices_ ? trkSelectors[i].getParameter<int32_t>("vtxNumber") : 0);
90  vertexCut_.push_back(useVertices_ ? trkSelectors[i].getParameter<std::string>("vertexCut") : nullptr);
91  // parameters for adapted optimal cuts on chi2 and primary vertex compatibility
92  res_par_.push_back(trkSelectors[i].getParameter<std::vector<double>>("res_par"));
93  chi2n_par_.push_back(trkSelectors[i].getParameter<double>("chi2n_par"));
94  chi2n_no1Dmod_par_.push_back(trkSelectors[i].getParameter<double>("chi2n_no1Dmod_par"));
95  d0_par1_.push_back(trkSelectors[i].getParameter<std::vector<double>>("d0_par1"));
96  dz_par1_.push_back(trkSelectors[i].getParameter<std::vector<double>>("dz_par1"));
97  d0_par2_.push_back(trkSelectors[i].getParameter<std::vector<double>>("d0_par2"));
98  dz_par2_.push_back(trkSelectors[i].getParameter<std::vector<double>>("dz_par2"));
99  // Boolean indicating if adapted primary vertex compatibility cuts are to be applied.
100  applyAdaptedPVCuts_.push_back(trkSelectors[i].getParameter<bool>("applyAdaptedPVCuts"));
101  // Impact parameter absolute cuts.
102  max_d0_.push_back(trkSelectors[i].getParameter<double>("max_d0"));
103  max_z0_.push_back(trkSelectors[i].getParameter<double>("max_z0"));
104  nSigmaZ_.push_back(trkSelectors[i].getParameter<double>("nSigmaZ"));
105  // Cuts on numbers of layers with hits/3D hits/lost hits.
106  min_layers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumberLayers"));
107  min_3Dlayers_.push_back(trkSelectors[i].getParameter<uint32_t>("minNumber3DLayers"));
108  max_lostLayers_.push_back(trkSelectors[i].getParameter<uint32_t>("maxNumberLostLayers"));
109  min_hits_bypass_.push_back(trkSelectors[i].getParameter<uint32_t>("minHitsToBypassChecks"));
110  // Flag to apply absolute cuts if no PV passes the selection
111  applyAbsCutsIfNoPV_.push_back(trkSelectors[i].getParameter<bool>("applyAbsCutsIfNoPV"));
112  keepAllTracks_.push_back(trkSelectors[i].getParameter<bool>("keepAllTracks"));
113  max_relpterr_.push_back(trkSelectors[i].getParameter<double>("max_relpterr"));
114  min_nhits_.push_back(trkSelectors[i].getParameter<uint32_t>("min_nhits"));
115  max_minMissHitOutOrIn_.push_back(trkSelectors[i].existsAs<int32_t>("max_minMissHitOutOrIn")
116  ? trkSelectors[i].getParameter<int32_t>("max_minMissHitOutOrIn")
117  : 99);
118  max_lostHitFraction_.push_back(trkSelectors[i].existsAs<double>("max_lostHitFraction")
119  ? trkSelectors[i].getParameter<double>("max_lostHitFraction")
120  : 1.0);
121  min_eta_.push_back(trkSelectors[i].existsAs<double>("min_eta") ? trkSelectors[i].getParameter<double>("min_eta")
122  : -9999);
123  max_eta_.push_back(trkSelectors[i].existsAs<double>("max_eta") ? trkSelectors[i].getParameter<double>("max_eta")
124  : 9999);
125 
126  setQualityBit_.push_back(false);
127  std::string qualityStr = trkSelectors[i].getParameter<std::string>("qualityBit");
128  if (!qualityStr.empty()) {
129  setQualityBit_[i] = true;
130  qualityToSet_[i] = TrackBase::qualityByName(trkSelectors[i].getParameter<std::string>("qualityBit"));
131  }
132 
133  if (setQualityBit_[i] && (qualityToSet_[i] == TrackBase::undefQuality))
134  throw cms::Exception("Configuration")
135  << "You can't set the quality bit " << trkSelectors[i].getParameter<std::string>("qualityBit")
136  << " as it is 'undefQuality' or unknown.\n";
137 
138  if (applyAbsCutsIfNoPV_[i]) {
139  max_d0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_d0NoPV"));
140  max_z0NoPV_.push_back(trkSelectors[i].getParameter<double>("max_z0NoPV"));
141  } else { //dummy values
142  max_d0NoPV_.push_back(0.);
143  max_z0NoPV_.push_back(0.);
144  }
145 
146  name_.push_back(trkSelectors[i].getParameter<std::string>("name"));
147 
148  preFilter_[i] = trkSelectors.size(); // no prefilter
149 
150  std::string pfName = trkSelectors[i].getParameter<std::string>("preFilterName");
151  if (!pfName.empty()) {
152  bool foundPF = false;
153  for (unsigned int j = 0; j < i; j++)
154  if (name_[j] == pfName) {
155  foundPF = true;
156  preFilter_[i] = j;
157  }
158  if (!foundPF)
159  throw cms::Exception("Configuration") << "Invalid prefilter name in MultiTrackSelector "
160  << trkSelectors[i].getParameter<std::string>("preFilterName");
161  }
162 
163  // produces<std::vector<int> >(name_[i]).setBranchAlias( name_[i] + "TrackQuals");
164  produces<edm::ValueMap<int>>(name_[i]).setBranchAlias(name_[i] + "TrackQuals");
165  produces<QualityMaskCollection>(name_[i]).setBranchAlias(name_[i] + "QualityMasks");
166  if (useAnyMVA_) {
167  bool thisMVA = false;
168  if (trkSelectors[i].exists("useMVA"))
169  thisMVA = trkSelectors[i].getParameter<bool>("useMVA");
170  useMVA_.push_back(thisMVA);
171  if (thisMVA) {
172  double minVal = -1;
173  if (trkSelectors[i].exists("minMVA"))
174  minVal = trkSelectors[i].getParameter<double>("minMVA");
175  min_MVA_.push_back(minVal);
176  mvaType_.push_back(trkSelectors[i].exists("mvaType") ? trkSelectors[i].getParameter<std::string>("mvaType")
177  : "Detached");
178  forestLabel_.push_back(trkSelectors[i].exists("GBRForestLabel")
179  ? trkSelectors[i].getParameter<std::string>("GBRForestLabel")
180  : "MVASelectorIter0");
181  useMVAonly_.push_back(trkSelectors[i].exists("useMVAonly") ? trkSelectors[i].getParameter<bool>("useMVAonly")
182  : false);
183  } else {
184  min_MVA_.push_back(-9999.0);
185  useMVAonly_.push_back(false);
186  mvaType_.push_back("Detached");
187  forestLabel_.push_back("MVASelectorIter0");
188  }
189  } else {
190  useMVA_.push_back(false);
191  useMVAonly_.push_back(false);
192  min_MVA_.push_back(-9999.0);
193  mvaType_.push_back("Detached");
194  forestLabel_.push_back("MVASelectorIter0");
195  }
196  }
197 
198  if (useForestFromDB_) {
199  forestToken_ = edm::vector_transform(forestLabel_, [this](auto const& label) {
200  return esConsumes<GBRForest, GBRWrapperRcd>(edm::ESInputTag("", label));
201  });
202  }
203 }
std::vector< uint32_t > min_nhits_
std::vector< double > max_d0NoPV_
std::vector< double > max_z0_
std::vector< uint32_t > min_layers_
Cuts on numbers of layers with hits/3D hits/lost hits.
std::vector< std::string > name_
std::vector< double > max_z0NoPV_
std::vector< bool > useMVA_
Quality qualityByName(std::string const &name)
std::vector< bool > keepAllTracks_
std::vector< uint32_t > min_hits_bypass_
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< bool > applyAdaptedPVCuts_
std::vector< double > chi2n_no1Dmod_par_
std::vector< double > nSigmaZ_
std::vector< std::string > mvaType_
std::vector< unsigned int > preFilter_
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
edm::EDGetTokenT< TrackingRecHitCollection > hSrc_
edm::EDGetTokenT< reco::VertexCollection > vertices_
std::vector< uint32_t > min_3Dlayers_
char const * label
std::vector< uint32_t > max_lostLayers_
std::vector< GBRForest * > forest_
std::vector< int32_t > max_lostHitFraction_
edm::EDGetTokenT< reco::BeamSpot > beamspot_
std::vector< double > min_eta_
std::vector< bool > setQualityBit_
do I have to set a quality bit?
std::vector< std::vector< double > > dz_par1_
std::vector< std::vector< double > > d0_par1_
std::vector< double > max_relpterr_
std::vector< int32_t > vtxNumber_
vertex cuts
std::vector< reco::TrackBase::TrackQuality > qualityToSet_
std::vector< edm::ESGetToken< GBRForest, GBRWrapperRcd > > forestToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< std::string > forestLabel_
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
std::vector< std::vector< double > > dz_par2_
std::vector< bool > applyAbsCutsIfNoPV_
std::vector< int32_t > max_minMissHitOutOrIn_
std::vector< double > min_MVA_
Log< level::Warning, false > LogWarning
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
std::vector< double > max_d0_
Impact parameter absolute cuts.
std::vector< bool > useMVAonly_
std::vector< std::vector< double > > res_par_
std::vector< double > max_eta_
std::vector< double > chi2n_par_
std::vector< std::vector< double > > d0_par2_
MultiTrackSelector::~MultiTrackSelector ( )
override

destructor

Definition at line 205 of file MultiTrackSelector.cc.

References forest_.

205  {
206  for (auto forest : forest_)
207  delete forest;
208 }
std::vector< GBRForest * > forest_

Member Function Documentation

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

Definition at line 210 of file MultiTrackSelector.cc.

References dbFileName_, forest_, forestLabel_, mps_fire::i, and useForestFromDB_.

210  {
211  if (!useForestFromDB_) {
212  TFile gbrfile(dbFileName_.c_str());
213  for (int i = 0; i < (int)forestLabel_.size(); i++) {
214  forest_[i] = (GBRForest*)gbrfile.Get(forestLabel_[i].c_str());
215  }
216  }
217 }
std::vector< GBRForest * > forest_
std::vector< std::string > forestLabel_
MultiTrackSelector::Point MultiTrackSelector::getBestVertex ( const reco::TrackBaseRef  track,
const reco::VertexCollection  vertices 
) const
protected

Definition at line 656 of file MultiTrackSelector.cc.

References PVValHelper::dz, reco::TrackBase::dz(), edm::RefToBase< T >::get(), AlCaHLTBitMon_ParallelJobs::p, and w().

Referenced by processMVA().

656  {
657  Point p(0, 0, -99999);
658  Point p_dz(0, 0, -99999);
659  float bestWeight = 0;
660  float dzmin = 10000;
661  bool weightMatch = false;
662  for (auto const& vertex : vertices) {
663  float w = vertex.trackWeight(track);
664  const Point& v_pos = vertex.position();
665  if (w > bestWeight) {
666  p = v_pos;
667  bestWeight = w;
668  weightMatch = true;
669  }
670  float dz = fabs(track.get()->dz(v_pos));
671  if (dz < dzmin) {
672  p_dz = v_pos;
673  dzmin = dz;
674  }
675  }
676  if (weightMatch)
677  return p;
678  else
679  return p_dz;
680 }
value_type const * get() const
Definition: RefToBase.h:209
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
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
T w() const
void MultiTrackSelector::processMVA ( edm::Event evt,
const edm::EventSetup es,
const reco::BeamSpot beamspot,
const reco::VertexCollection vertices,
int  selIndex,
std::vector< float > &  mvaVals_,
bool  writeIt = false 
) const
protected

Definition at line 532 of file MultiTrackSelector.cc.

References cms::cuda::assert(), reco::TrackBase::chi2(), HLT_FULL_cff::chi2, chi2n, TrackingRecHit::dimension(), reco::TrackBase::dxy(), reco::TrackBase::dz(), reco::TrackBase::eta(), reco::Track::extra(), validate-o2o-wbm::f, forest_, forestToken_, getBestVertex(), edm::Event::getByToken(), GBRForest::GetClassifier(), edm::EventSetup::getData(), reco::TrackBase::hitPattern(), hSrc_, SiStripPI::max, SiStripPI::min, reco::HitPattern::MISSING_INNER_HITS, reco::HitPattern::MISSING_OUTER_HITS, eostools::move(), mvaType_, ndof, reco::TrackBase::ndof(), reco::TrackBase::normalizedChi2(), reco::TrackBase::numberOfLostHits(), reco::HitPattern::numberOfLostTrackerHits(), reco::TrackBase::numberOfValidHits(), reco::HitPattern::numberOfValidStripLayersWithMonoAndStereo(), reco::HitPattern::pixelLayersWithMeasurement(), reco::BeamSpot::position(), reco::TrackBase::pt(), reco::TrackBase::ptError(), edm::Event::put(), reco::Track::recHitsSize(), dt_dqm_sourceclient_common_cff::reco, src_, HiCentrality_cfi::srcTracks, reco::HitPattern::TRACK_HITS, reco::HitPattern::trackerLayersWithMeasurement(), reco::HitPattern::trackerLayersWithoutMeasurement(), useAnyMVA_, useForestFromDB_, and useMVA_.

Referenced by AnalyticalTrackSelector::run(), and run().

538  {
539  using namespace std;
540  using namespace edm;
541  using namespace reco;
542 
543  // Get tracks
544  Handle<TrackCollection> hSrcTrack;
545  evt.getByToken(src_, hSrcTrack);
546  const TrackCollection& srcTracks(*hSrcTrack);
547  RefToBaseProd<Track> rtbpTrackCollection(hSrcTrack);
548  assert(mvaVals_.size() == srcTracks.size());
549 
550  // get hits in track..
552  evt.getByToken(hSrc_, hSrcHits);
553  const TrackingRecHitCollection& srcHits(*hSrcHits);
554 
555  auto mvaValValueMap = std::make_unique<edm::ValueMap<float>>();
556  edm::ValueMap<float>::Filler mvaFiller(*mvaValValueMap);
557 
558  if (!useAnyMVA_ && writeIt) {
559  // mvaVals_ already initalized...
560  mvaFiller.insert(hSrcTrack, mvaVals_.begin(), mvaVals_.end());
561  mvaFiller.fill();
562  evt.put(std::move(mvaValValueMap), "MVAVals");
563  auto mvas = std::make_unique<MVACollection>(mvaVals_.begin(), mvaVals_.end());
564  evt.put(std::move(mvas), "MVAValues");
565  return;
566  }
567 
568  if (!useMVA_[selIndex] && !writeIt)
569  return;
570 
571  size_t current = 0;
572  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
573  const Track& trk = *it;
574  RefToBase<Track> trackRef(rtbpTrackCollection, current);
575  auto tmva_ndof_ = trk.ndof();
576  auto tmva_nlayers_ = trk.hitPattern().trackerLayersWithMeasurement();
577  auto tmva_nlayers3D_ =
580  float chi2n = trk.normalizedChi2();
581  float chi2n_no1Dmod = chi2n;
582 
583  int count1dhits = 0;
584  auto ith = trk.extra()->firstRecHit();
585  auto edh = ith + trk.recHitsSize();
586  for (; ith < edh; ++ith) {
587  const TrackingRecHit& hit = srcHits[ith];
588  if (hit.dimension() == 1)
589  ++count1dhits;
590  }
591  if (count1dhits > 0) {
592  float chi2 = trk.chi2();
593  float ndof = trk.ndof();
594  chi2n = (chi2 + count1dhits) / float(ndof + count1dhits);
595  }
596  auto tmva_chi2n_ = chi2n;
597  auto tmva_chi2n_no1dmod_ = chi2n_no1Dmod;
598  auto tmva_eta_ = trk.eta();
599  auto tmva_relpterr_ = float(trk.ptError()) / std::max(float(trk.pt()), 0.000001f);
600  auto tmva_nhits_ = trk.numberOfValidHits();
603  auto tmva_minlost_ = std::min(lostIn, lostOut);
604  auto tmva_lostmidfrac_ = trk.numberOfLostHits() / (trk.numberOfValidHits() + trk.numberOfLostHits());
605  auto tmva_absd0_ = fabs(-trk.dxy(beamspot.position()));
606  auto tmva_absdz_ = fabs(trk.dz(beamspot.position()));
607  Point bestVertex = getBestVertex(trackRef, vertices);
608  auto tmva_absd0PV_ = fabs(trk.dxy(bestVertex));
609  auto tmva_absdzPV_ = fabs(trk.dz(bestVertex));
610  auto tmva_pt_ = trk.pt();
611 
612  GBRForest const* forest = forest_[selIndex];
613  if (useForestFromDB_) {
614  forest = &es.getData(forestToken_[selIndex]);
615  }
616 
617  float gbrVals_[16];
618  gbrVals_[0] = tmva_pt_;
619  gbrVals_[1] = tmva_lostmidfrac_;
620  gbrVals_[2] = tmva_minlost_;
621  gbrVals_[3] = tmva_nhits_;
622  gbrVals_[4] = tmva_relpterr_;
623  gbrVals_[5] = tmva_eta_;
624  gbrVals_[6] = tmva_chi2n_no1dmod_;
625  gbrVals_[7] = tmva_chi2n_;
626  gbrVals_[8] = tmva_nlayerslost_;
627  gbrVals_[9] = tmva_nlayers3D_;
628  gbrVals_[10] = tmva_nlayers_;
629  gbrVals_[11] = tmva_ndof_;
630  gbrVals_[12] = tmva_absd0PV_;
631  gbrVals_[13] = tmva_absdzPV_;
632  gbrVals_[14] = tmva_absdz_;
633  gbrVals_[15] = tmva_absd0_;
634 
635  if (mvaType_[selIndex] == "Prompt") {
636  auto gbrVal = forest->GetClassifier(gbrVals_);
637  mvaVals_[current] = gbrVal;
638  } else {
639  float detachedGbrVals_[12];
640  for (int jjj = 0; jjj < 12; jjj++)
641  detachedGbrVals_[jjj] = gbrVals_[jjj];
642  auto gbrVal = forest->GetClassifier(detachedGbrVals_);
643  mvaVals_[current] = gbrVal;
644  }
645  }
646 
647  if (writeIt) {
648  mvaFiller.insert(hSrcTrack, mvaVals_.begin(), mvaVals_.end());
649  mvaFiller.fill();
650  evt.put(std::move(mvaValValueMap), "MVAVals");
651  auto mvas = std::make_unique<MVACollection>(mvaVals_.begin(), mvaVals_.end());
652  evt.put(std::move(mvas), "MVAValues");
653  }
654 }
virtual int dimension() const =0
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:139
std::vector< bool > useMVA_
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
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
Point getBestVertex(const reco::TrackBaseRef, const reco::VertexCollection) const
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:801
std::vector< std::string > mvaType_
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:513
edm::EDGetTokenT< TrackingRecHitCollection > hSrc_
assert(be >=bs)
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:534
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:893
bool getData(T &iHolder) const
Definition: EventSetup.h:122
std::vector< GBRForest * > forest_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:369
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:590
double pt() const
track transverse momentum
Definition: TrackBase.h:637
def move
Definition: eostools.py:511
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
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
const HitPattern & hitPattern() const
Access the hit pattern, indicating in which Tracker layers the track has hits.
Definition: TrackBase.h:504
std::vector< edm::ESGetToken< GBRForest, GBRWrapperRcd > > forestToken_
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:553
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
const Point & position() const
position
Definition: BeamSpot.h:59
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
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
double GetClassifier(const float *vector) const
Definition: GBRForest.h:33
void MultiTrackSelector::produce ( edm::Event evt,
const edm::EventSetup es 
)
inlinefinalprotected

process one event

Definition at line 59 of file MultiTrackSelector.h.

References submitPVValidationJobs::run.

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

Reimplemented in AnalyticalTrackSelector.

Definition at line 219 of file MultiTrackSelector.cc.

References beamspot_, edm::Event::getByToken(), hSrc_, mps_fire::i, dqmiolumiharvest::j, keepAllTracks_, LogTrace, loose, SiStripPI::max, eostools::move(), name_, convertSQLiteXML::ok, preFilter_, processMVA(), reco::TrackBase::pt(), edm::Event::put(), submitPVResolutionJobs::q, reco::TrackBase::qualityMask(), qualityToSet_, dt_dqm_sourceclient_common_cff::reco, select(), selectVertices(), setQualityBit_, src_, HiCentrality_cfi::srcTracks, useAnyMVA_, useVertices_, and vertices_.

219  {
220  using namespace std;
221  using namespace edm;
222  using namespace reco;
223 
224  // Get tracks
225  Handle<TrackCollection> hSrcTrack;
226  evt.getByToken(src_, hSrcTrack);
227 
228  const TrackCollection& srcTracks(*hSrcTrack);
229  if (hSrcTrack.failedToGet())
230  edm::LogWarning("MultiTrackSelector") << "could not get Track collection";
231 
232  // get hits in track..
234  evt.getByToken(hSrc_, hSrcHits);
235  const TrackingRecHitCollection& srcHits(*hSrcHits);
236 
237  // looking for the beam spot
239  evt.getByToken(beamspot_, hBsp);
240  const reco::BeamSpot& vertexBeamSpot(*hBsp);
241 
242  // Select good primary vertices for use in subsequent track selection
244  if (useVertices_) {
245  evt.getByToken(vertices_, hVtx);
246  if (hVtx.failedToGet())
247  edm::LogWarning("MultiTrackSelector") << "could not get Vertex collection";
248  }
249 
250  unsigned int trkSize = srcTracks.size();
251  std::vector<int> selTracksSave(qualityToSet_.size() * trkSize, 0);
252 
253  std::vector<Point> points;
254  std::vector<float> vterr, vzerr;
255  if (useVertices_)
256  selectVertices(0, *hVtx, points, vterr, vzerr);
257  //auto vtxP = points.empty() ? vertexBeamSpot.position() : points[0]; // rare, very rare, still happens!
258  for (unsigned int i = 0; i < qualityToSet_.size(); i++) {
259  std::vector<float> mvaVals_(srcTracks.size(), -99.f);
260  processMVA(evt, es, vertexBeamSpot, *(hVtx.product()), i, mvaVals_, i == 0 ? true : false);
261  std::vector<int> selTracks(trkSize, 0);
262  auto selTracksValueMap = std::make_unique<edm::ValueMap<int>>();
263  edm::ValueMap<int>::Filler filler(*selTracksValueMap);
264 
265  if (useVertices_)
266  selectVertices(i, *hVtx, points, vterr, vzerr);
267 
268  // Loop over tracks
269  size_t current = 0;
270  for (TrackCollection::const_iterator it = srcTracks.begin(), ed = srcTracks.end(); it != ed; ++it, ++current) {
271  const Track& trk = *it;
272  // Check if this track passes cuts
273 
274  LogTrace("TrackSelection") << "ready to check track with pt=" << trk.pt();
275 
276  //already removed
277  bool ok = true;
278  if (preFilter_[i] < i && selTracksSave[preFilter_[i] * trkSize + current] < 0) {
279  selTracks[current] = -1;
280  ok = false;
281  if (!keepAllTracks_[i])
282  continue;
283  } else {
284  float mvaVal = 0;
285  if (useAnyMVA_)
286  mvaVal = mvaVals_[current];
287  ok = select(i, vertexBeamSpot, srcHits, trk, points, vterr, vzerr, mvaVal);
288  if (!ok) {
289  LogTrace("TrackSelection") << "track with pt=" << trk.pt() << " NOT selected";
290  if (!keepAllTracks_[i]) {
291  selTracks[current] = -1;
292  continue;
293  }
294  } else
295  LogTrace("TrackSelection") << "track with pt=" << trk.pt() << " selected";
296  }
297 
298  if (preFilter_[i] < i) {
299  selTracks[current] = selTracksSave[preFilter_[i] * trkSize + current];
300  } else {
301  selTracks[current] = trk.qualityMask();
302  }
303  if (ok && setQualityBit_[i]) {
304  selTracks[current] = (selTracks[current] | (1 << qualityToSet_[i]));
305  if (qualityToSet_[i] == TrackBase::tight) {
306  selTracks[current] = (selTracks[current] | (1 << TrackBase::loose));
307  } else if (qualityToSet_[i] == TrackBase::highPurity) {
308  selTracks[current] = (selTracks[current] | (1 << TrackBase::loose));
309  selTracks[current] = (selTracks[current] | (1 << TrackBase::tight));
310  }
311  if (!points.empty()) {
312  if (qualityToSet_[i] == TrackBase::loose) {
313  selTracks[current] = (selTracks[current] | (1 << TrackBase::looseSetWithPV));
314  } else if (qualityToSet_[i] == TrackBase::highPurity) {
315  selTracks[current] = (selTracks[current] | (1 << TrackBase::looseSetWithPV));
316  selTracks[current] = (selTracks[current] | (1 << TrackBase::highPuritySetWithPV));
317  }
318  }
319  }
320  }
321  for (unsigned int j = 0; j < trkSize; j++)
322  selTracksSave[j + i * trkSize] = selTracks[j];
323  filler.insert(hSrcTrack, selTracks.begin(), selTracks.end());
324  filler.fill();
325 
326  // evt.put(std::move(selTracks),name_[i]);
327  evt.put(std::move(selTracksValueMap), name_[i]);
328  for (auto& q : selTracks)
329  q = std::max(q, 0);
330  auto quals = std::make_unique<QualityMaskCollection>(selTracks.begin(), selTracks.end());
331  evt.put(std::move(quals), name_[i]);
332  }
333 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
std::vector< std::string > name_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
std::vector< bool > keepAllTracks_
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
std::vector< Track > TrackCollection
collection of Tracks
Definition: TrackFwd.h:14
std::vector< unsigned int > preFilter_
edm::EDGetTokenT< TrackingRecHitCollection > hSrc_
edm::EDGetTokenT< reco::VertexCollection > vertices_
#define LogTrace(id)
edm::EDGetTokenT< reco::BeamSpot > beamspot_
std::vector< bool > setQualityBit_
do I have to set a quality bit?
double pt() const
track transverse momentum
Definition: TrackBase.h:637
def move
Definition: eostools.py:511
int qualityMask() const
Definition: TrackBase.h:843
void selectVertices(unsigned int tsNum, const reco::VertexCollection &vtxs, std::vector< Point > &points, std::vector< float > &vterr, std::vector< float > &vzerr) const
std::vector< reco::TrackBase::TrackQuality > qualityToSet_
constexpr auto loose
edm::EDGetTokenT< reco::TrackCollection > src_
source collection label
void processMVA(edm::Event &evt, const edm::EventSetup &es, const reco::BeamSpot &beamspot, const reco::VertexCollection &vertices, int selIndex, std::vector< float > &mvaVals_, bool writeIt=false) const
bool MultiTrackSelector::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 335 of file MultiTrackSelector.cc.

References funct::abs(), applyAbsCutsIfNoPV_, applyAdaptedPVCuts_, reco::TrackBase::chi2(), HLT_FULL_cff::chi2, chi2n, chi2n_no1Dmod_par_, chi2n_par_, d0, d0_par1_, d0_par2_, reco::TrackBase::d0Error(), TrackingRecHit::dimension(), reco::TrackBase::dxy(), PVValHelper::dz, reco::TrackBase::dz(), dz_par1_, dz_par2_, reco::TrackBase::dzError(), dataset::end, PVValHelper::eta, reco::TrackBase::eta(), reco::Track::extra(), validate-o2o-wbm::f, 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(), reco::HitPattern::pixelLayersWithMeasurement(), point, reco::BeamSpot::position(), funct::pow(), DiDispStaMuonMonitor_cfi::pt, reco::TrackBase::pt(), reco::TrackBase::ptError(), reco::Track::recHitsSize(), relpterr, res_par_, reco::BeamSpot::sigmaZ(), mathSSE::sqrt(), reco::HitPattern::TRACK_HITS, reco::HitPattern::trackerLayersWithMeasurement(), reco::HitPattern::trackerLayersWithoutMeasurement(), useAnyMVA_, useMVA_, useMVAonly_, and useVtxError_.

Referenced by AnalyticalTrackSelector::run(), and run().

342  {
343  // Decide if the given track passes selection cuts.
344 
345  using namespace std;
346 
347  //cuts on number of valid hits
348  auto nhits = tk.numberOfValidHits();
349  if (nhits >= min_hits_bypass_[tsNum])
350  return true;
351  if (nhits < min_nhits_[tsNum])
352  return false;
353 
354  if (tk.ndof() < 1E-5)
355  return false;
356 
358  //Adding the MVA selection before any other cut//
360  if (useAnyMVA_ && useMVA_[tsNum]) {
361  if (useMVAonly_[tsNum])
362  return mvaVal > min_MVA_[tsNum];
363  if (mvaVal < min_MVA_[tsNum])
364  return false;
365  }
367  //End of MVA selection section//
369 
370  // Cuts on numbers of layers with hits/3D hits/lost hits.
372  uint32_t nlayers3D =
375  LogDebug("TrackSelection") << "cuts on nlayers: " << nlayers << " " << nlayers3D << " " << nlayersLost << " vs "
376  << min_layers_[tsNum] << " " << min_3Dlayers_[tsNum] << " " << max_lostLayers_[tsNum];
377  if (nlayers < min_layers_[tsNum])
378  return false;
379  if (nlayers3D < min_3Dlayers_[tsNum])
380  return false;
381  if (nlayersLost > max_lostLayers_[tsNum])
382  return false;
383  LogTrace("TrackSelection") << "cuts on nlayers passed";
384 
385  float chi2n = tk.normalizedChi2();
386  float chi2n_no1Dmod = chi2n;
387 
388  int count1dhits = 0;
389  auto ith = tk.extra()->firstRecHit();
390  auto edh = ith + tk.recHitsSize();
391  for (; ith < edh; ++ith) {
392  const TrackingRecHit& hit = recHits[ith];
393  if (hit.dimension() == 1)
394  ++count1dhits;
395  }
396  if (count1dhits > 0) {
397  float chi2 = tk.chi2();
398  float ndof = tk.ndof();
399  chi2n = (chi2 + count1dhits) / float(ndof + count1dhits);
400  }
401  // For each 1D rechit, the chi^2 and ndof is increased by one. This is a way of retaining approximately
402  // the same normalized chi^2 distribution as with 2D rechits.
403  if (chi2n > chi2n_par_[tsNum] * nlayers)
404  return false;
405 
406  if (chi2n_no1Dmod > chi2n_no1Dmod_par_[tsNum] * nlayers)
407  return false;
408 
409  // Get track parameters
410  float pt = std::max(float(tk.pt()), 0.000001f);
411  float eta = tk.eta();
412  if (eta < min_eta_[tsNum] || eta > max_eta_[tsNum])
413  return false;
414 
415  //cuts on relative error on pt
416  float relpterr = float(tk.ptError()) / pt;
417  if (relpterr > max_relpterr_[tsNum])
418  return false;
419 
422  int minLost = std::min(lostIn, lostOut);
423  if (minLost > max_minMissHitOutOrIn_[tsNum])
424  return false;
425  float lostMidFrac = tk.numberOfLostHits() / (tk.numberOfValidHits() + tk.numberOfLostHits());
426  if (lostMidFrac > max_lostHitFraction_[tsNum])
427  return false;
428 
429  //other track parameters
430  float d0 = -tk.dxy(vertexBeamSpot.position()), d0E = tk.d0Error(), dz = tk.dz(vertexBeamSpot.position()),
431  dzE = tk.dzError();
432 
433  // parametrized d0 resolution for the track pt
434  float nomd0E = sqrt(res_par_[tsNum][0] * res_par_[tsNum][0] + (res_par_[tsNum][1] / pt) * (res_par_[tsNum][1] / pt));
435  // parametrized z0 resolution for the track pt and eta
436  float nomdzE = nomd0E * (std::cosh(eta));
437 
438  float dzCut = std::min(powN(dz_par1_[tsNum][0] * nlayers, int(dz_par1_[tsNum][1] + 0.5)) * nomdzE,
439  powN(dz_par2_[tsNum][0] * nlayers, int(dz_par2_[tsNum][1] + 0.5)) * dzE);
440  float d0Cut = std::min(powN(d0_par1_[tsNum][0] * nlayers, int(d0_par1_[tsNum][1] + 0.5)) * nomd0E,
441  powN(d0_par2_[tsNum][0] * nlayers, int(d0_par2_[tsNum][1] + 0.5)) * d0E);
442 
443  // ---- PrimaryVertex compatibility cut
444  bool primaryVertexZCompatibility(false);
445  bool primaryVertexD0Compatibility(false);
446 
447  if (points.empty()) { //If not primaryVertices are reconstructed, check just the compatibility with the BS
448  //z0 within (n sigma + dzCut) of the beam spot z, if no good vertex is found
449  if (abs(dz) < hypot(vertexBeamSpot.sigmaZ() * nSigmaZ_[tsNum], dzCut))
450  primaryVertexZCompatibility = true;
451  // d0 compatibility with beam line
452  if (abs(d0) < d0Cut)
453  primaryVertexD0Compatibility = true;
454  }
455 
456  int iv = 0;
457  for (std::vector<Point>::const_iterator point = points.begin(), end = points.end(); point != end; ++point) {
458  LogTrace("TrackSelection") << "Test track w.r.t. vertex with z position " << point->z();
459  if (primaryVertexZCompatibility && primaryVertexD0Compatibility)
460  break;
461  float dzPV = tk.dz(*point); //re-evaluate the dz with respect to the vertex position
462  float d0PV = tk.dxy(*point); //re-evaluate the dxy with respect to the vertex position
463  if (useVtxError_) {
464  float dzErrPV = std::sqrt(dzE * dzE + vzerr[iv] * vzerr[iv]); // include vertex error in z
465  float d0ErrPV = std::sqrt(d0E * d0E + vterr[iv] * vterr[iv]); // include vertex error in xy
466  iv++;
467  if (abs(dzPV) < dz_par1_[tsNum][0] * pow(nlayers, dz_par1_[tsNum][1]) * nomdzE &&
468  abs(dzPV) < dz_par2_[tsNum][0] * pow(nlayers, dz_par2_[tsNum][1]) * dzErrPV && abs(dzPV) < max_z0_[tsNum])
469  primaryVertexZCompatibility = true;
470  if (abs(d0PV) < d0_par1_[tsNum][0] * pow(nlayers, d0_par1_[tsNum][1]) * nomd0E &&
471  abs(d0PV) < d0_par2_[tsNum][0] * pow(nlayers, d0_par2_[tsNum][1]) * d0ErrPV && abs(d0PV) < max_d0_[tsNum])
472  primaryVertexD0Compatibility = true;
473  } else {
474  if (abs(dzPV) < dzCut)
475  primaryVertexZCompatibility = true;
476  if (abs(d0PV) < d0Cut)
477  primaryVertexD0Compatibility = true;
478  }
479  LogTrace("TrackSelection") << "distances " << dzPV << " " << d0PV << " vs " << dzCut << " " << d0Cut;
480  }
481 
482  if (points.empty() && applyAbsCutsIfNoPV_[tsNum]) {
483  if (abs(dz) > max_z0NoPV_[tsNum] || abs(d0) > max_d0NoPV_[tsNum])
484  return false;
485  } else {
486  // Absolute cuts on all tracks impact parameters with respect to beam-spot.
487  // If BS is not compatible, verify if at least the reco-vertex is compatible (useful for incorrect BS settings)
488  if (abs(d0) > max_d0_[tsNum] && !primaryVertexD0Compatibility)
489  return false;
490  LogTrace("TrackSelection") << "absolute cuts on d0 passed";
491  if (abs(dz) > max_z0_[tsNum] && !primaryVertexZCompatibility)
492  return false;
493  LogTrace("TrackSelection") << "absolute cuts on dz passed";
494  }
495 
496  LogTrace("TrackSelection") << "cuts on PV: apply adapted PV cuts? " << applyAdaptedPVCuts_[tsNum]
497  << " d0 compatibility? " << primaryVertexD0Compatibility << " z compatibility? "
498  << primaryVertexZCompatibility;
499 
500  if (applyAdaptedPVCuts_[tsNum]) {
501  return (primaryVertexD0Compatibility && primaryVertexZCompatibility);
502  } else {
503  return true;
504  }
505 }
std::vector< uint32_t > min_nhits_
virtual int dimension() const =0
std::vector< double > max_d0NoPV_
double d0Error() const
error on d0
Definition: TrackBase.h:772
std::vector< double > max_z0_
std::vector< uint32_t > min_layers_
Cuts on numbers of layers with hits/3D hits/lost hits.
std::vector< double > max_z0NoPV_
int32_t *__restrict__ iv
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:139
std::vector< bool > useMVA_
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
size_t recHitsSize() const
Get number of RecHits. (Warning, this includes invalid hits, which are not physical hits)...
Definition: Track.h:97
std::vector< uint32_t > min_hits_bypass_
std::vector< bool > applyAdaptedPVCuts_
unsigned short numberOfLostHits() const
number of cases where track crossed a layer without getting a hit.
Definition: TrackBase.h:801
std::vector< double > chi2n_no1Dmod_par_
std::vector< double > nSigmaZ_
int pixelLayersWithMeasurement() const
Definition: HitPattern.cc:513
int trackerLayersWithMeasurement() const
Definition: HitPattern.cc:534
#define LogTrace(id)
std::vector< uint32_t > min_3Dlayers_
int numberOfLostTrackerHits(HitCategory category) const
Definition: HitPattern.h:893
std::vector< uint32_t > max_lostLayers_
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
int numberOfValidStripLayersWithMonoAndStereo(uint16_t stripdet, uint16_t layer) const
Definition: HitPattern.cc:369
std::vector< int32_t > max_lostHitFraction_
std::vector< double > min_eta_
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:590
T sqrt(T t)
Definition: SSEVec.h:19
double pt() const
track transverse momentum
Definition: TrackBase.h:637
std::vector< std::vector< double > > dz_par1_
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:754
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
unsigned short numberOfValidHits() const
number of valid hits found
Definition: TrackBase.h:798
std::vector< std::vector< double > > d0_par1_
std::vector< double > max_relpterr_
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
static constexpr float d0
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
int trackerLayersWithoutMeasurement(HitCategory category) const
Definition: HitPattern.cc:553
string end
Definition: dataset.py:937
std::vector< std::vector< double > > dz_par2_
std::vector< bool > applyAbsCutsIfNoPV_
const Point & position() const
position
Definition: BeamSpot.h:59
std::vector< int32_t > max_minMissHitOutOrIn_
std::vector< double > min_MVA_
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
std::vector< double > max_d0_
Impact parameter absolute cuts.
std::vector< bool > useMVAonly_
std::vector< std::vector< double > > res_par_
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
std::vector< double > max_eta_
std::vector< double > chi2n_par_
#define LogDebug(id)
std::vector< std::vector< double > > d0_par2_
void MultiTrackSelector::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 507 of file MultiTrackSelector.cc.

References LogDebug, LogTrace, dt_dqm_sourceclient_common_cff::reco, mathSSE::sqrt(), vertexCut_, and vtxNumber_.

Referenced by AnalyticalTrackSelector::run(), and run().

511  {
512  // Select good primary vertices
513  using namespace reco;
514  int32_t toTake = vtxNumber_[tsNum];
515  for (VertexCollection::const_iterator it = vtxs.begin(), ed = vtxs.end(); it != ed; ++it) {
516  LogDebug("SelectVertex") << " select vertex with z position " << it->z() << " " << it->chi2() << " " << it->ndof()
517  << " " << TMath::Prob(it->chi2(), static_cast<int32_t>(it->ndof()));
518  Vertex vtx = *it;
519  bool pass = vertexCut_[tsNum](vtx);
520  if (pass) {
521  points.push_back(it->position());
522  vterr.push_back(sqrt(it->yError() * it->xError()));
523  vzerr.push_back(it->zError());
524  LogTrace("SelectVertex") << " SELECTED vertex with z position " << it->z();
525  toTake--;
526  if (toTake == 0)
527  break;
528  }
529  }
530 }
#define LogTrace(id)
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< int32_t > vtxNumber_
vertex cuts
std::vector< StringCutObjectSelector< reco::Vertex > > vertexCut_
#define LogDebug(id)

Member Data Documentation

std::vector<bool> MultiTrackSelector::applyAbsCutsIfNoPV_
protected
std::vector<bool> MultiTrackSelector::applyAdaptedPVCuts_
protected
edm::EDGetTokenT<reco::BeamSpot> MultiTrackSelector::beamspot_
protected
std::vector<double> MultiTrackSelector::chi2n_no1Dmod_par_
protected
std::vector<double> MultiTrackSelector::chi2n_par_
protected
std::vector<std::vector<double> > MultiTrackSelector::d0_par1_
protected
std::vector<std::vector<double> > MultiTrackSelector::d0_par2_
protected
std::string MultiTrackSelector::dbFileName_
protected

Definition at line 158 of file MultiTrackSelector.h.

Referenced by beginStream(), and MultiTrackSelector().

std::vector<std::vector<double> > MultiTrackSelector::dz_par1_
protected
std::vector<std::vector<double> > MultiTrackSelector::dz_par2_
protected
std::vector<GBRForest *> MultiTrackSelector::forest_
protected
std::vector<std::string> MultiTrackSelector::forestLabel_
protected
std::vector<edm::ESGetToken<GBRForest, GBRWrapperRcd> > MultiTrackSelector::forestToken_
protected

Definition at line 155 of file MultiTrackSelector.h.

Referenced by MultiTrackSelector(), and processMVA().

edm::EDGetTokenT<TrackingRecHitCollection> MultiTrackSelector::hSrc_
protected
std::vector<bool> MultiTrackSelector::keepAllTracks_
protected
std::vector<double> MultiTrackSelector::max_d0_
protected

Impact parameter absolute cuts.

Definition at line 116 of file MultiTrackSelector.h.

Referenced by AnalyticalTrackSelector::AnalyticalTrackSelector(), MultiTrackSelector(), and select().

std::vector<double> MultiTrackSelector::max_d0NoPV_
protected
std::vector<double> MultiTrackSelector::max_eta_
protected
std::vector<int32_t> MultiTrackSelector::max_lostHitFraction_
protected
std::vector<uint32_t> MultiTrackSelector::max_lostLayers_
protected
std::vector<int32_t> MultiTrackSelector::max_minMissHitOutOrIn_
protected
std::vector<double> MultiTrackSelector::max_relpterr_
protected
std::vector<double> MultiTrackSelector::max_z0_
protected
std::vector<double> MultiTrackSelector::max_z0NoPV_
protected
std::vector<uint32_t> MultiTrackSelector::min_3Dlayers_
protected
std::vector<double> MultiTrackSelector::min_eta_
protected
std::vector<uint32_t> MultiTrackSelector::min_hits_bypass_
protected
std::vector<uint32_t> MultiTrackSelector::min_layers_
protected

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

Definition at line 121 of file MultiTrackSelector.h.

Referenced by AnalyticalTrackSelector::AnalyticalTrackSelector(), MultiTrackSelector(), and select().

std::vector<double> MultiTrackSelector::min_MVA_
protected
std::vector<uint32_t> MultiTrackSelector::min_nhits_
protected
std::vector<std::string> MultiTrackSelector::mvaType_
protected
std::vector<std::string> MultiTrackSelector::name_
protected

Definition at line 145 of file MultiTrackSelector.h.

Referenced by MultiTrackSelector(), and run().

std::vector<double> MultiTrackSelector::nSigmaZ_
protected
std::vector<unsigned int> MultiTrackSelector::preFilter_
protected
std::vector<reco::TrackBase::TrackQuality> MultiTrackSelector::qualityToSet_
protected
std::vector<std::vector<double> > MultiTrackSelector::res_par_
protected
std::vector<bool> MultiTrackSelector::setQualityBit_
protected

do I have to set a quality bit?

Definition at line 96 of file MultiTrackSelector.h.

Referenced by AnalyticalTrackSelector::AnalyticalTrackSelector(), MultiTrackSelector(), AnalyticalTrackSelector::run(), and run().

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

source collection label

Definition at line 87 of file MultiTrackSelector.h.

Referenced by AnalyticalTrackSelector::AnalyticalTrackSelector(), processMVA(), AnalyticalTrackSelector::run(), and run().

bool MultiTrackSelector::useAnyMVA_
protected
bool MultiTrackSelector::useForestFromDB_
protected

Definition at line 157 of file MultiTrackSelector.h.

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

std::vector<bool> MultiTrackSelector::useMVA_
protected
std::vector<bool> MultiTrackSelector::useMVAonly_
protected
bool MultiTrackSelector::useVertices_
protected
bool MultiTrackSelector::useVtxError_
protected
std::vector<StringCutObjectSelector<reco::Vertex> > MultiTrackSelector::vertexCut_
protected
edm::EDGetTokenT<reco::VertexCollection> MultiTrackSelector::vertices_
protected
std::vector<int32_t> MultiTrackSelector::vtxNumber_
protected