CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
SeedMvaEstimator Class Reference

#include <SeedMvaEstimator.h>

Public Member Functions

double computeMva (const TrajectorySeed &, const GlobalVector &, const l1t::MuonBxCollection &, const reco::RecoChargedCandidateCollection &) const
 
 SeedMvaEstimator (const edm::FileInPath &weightsfile, const std::vector< double > &scale_mean, const std::vector< double > &scale_std, const bool isFromL1, const int minL1Qual)
 
 ~SeedMvaEstimator ()
 

Private Member Functions

void getL1MuonVariables (const GlobalVector &, const l1t::MuonBxCollection &, float &, float &) const
 
void getL2MuonVariables (const GlobalVector &, const reco::RecoChargedCandidateCollection &, float &, float &) const
 

Private Attributes

std::unique_ptr< const GBRForestgbrForest_
 
const bool isFromL1_
 
const int minL1Qual_
 
const std::vector< double > scale_mean_
 
const std::vector< double > scale_std_
 

Detailed Description

Definition at line 23 of file SeedMvaEstimator.h.

Constructor & Destructor Documentation

◆ SeedMvaEstimator()

SeedMvaEstimator::SeedMvaEstimator ( const edm::FileInPath weightsfile,
const std::vector< double > &  scale_mean,
const std::vector< double > &  scale_std,
const bool  isFromL1,
const int  minL1Qual 
)

Definition at line 9 of file SeedMvaEstimator.cc.

References createGBRForest(), and gbrForest_.

14  : scale_mean_(scale_mean), scale_std_(scale_std), isFromL1_(isFromL1), minL1Qual_(minL1Qual) {
15  gbrForest_ = createGBRForest(weightsfile);
16 }
const bool isFromL1_
std::unique_ptr< const GBRForest > gbrForest_
const std::vector< double > scale_mean_
const std::vector< double > scale_std_
std::unique_ptr< const GBRForest > createGBRForest(const std::string &weightsFile)

◆ ~SeedMvaEstimator()

SeedMvaEstimator::~SeedMvaEstimator ( )

Definition at line 18 of file SeedMvaEstimator.cc.

18 {}

Member Function Documentation

◆ computeMva()

double SeedMvaEstimator::computeMva ( const TrajectorySeed seed,
const GlobalVector global_p,
const l1t::MuonBxCollection l1Muons,
const reco::RecoChargedCandidateCollection l2Muons 
) const

Definition at line 72 of file SeedMvaEstimator.cc.

References ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), gbrForest_, getL1MuonVariables(), getL2MuonVariables(), isFromL1_, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::iv, Phase2::kLast, Phase2::kTsosDxdz, Phase2::kTsosDydz, Phase2::kTsosErr0, Phase2::kTsosErr2, Phase2::kTsosErr5, Phase2::kTsosQbp, scale_mean_, scale_std_, fileCollector::seed, mathSSE::sqrt(), and trigObjTnPSource_cfi::var.

75  {
76  static constexpr float initDRdPhi(99999.);
77  auto kLast = isFromL1_ ? kLastL1 : kLastL2;
78  float var[kLast];
79 
80  var[kTsosErr0] = seed.startingState().error(0);
81  var[kTsosErr2] = seed.startingState().error(2);
82  var[kTsosErr5] = seed.startingState().error(5);
83  var[kTsosDxdz] = seed.startingState().parameters().dxdz();
84  var[kTsosDydz] = seed.startingState().parameters().dydz();
85  var[kTsosQbp] = seed.startingState().parameters().qbp();
86 
87  float dR2dRL1SeedP = initDRdPhi;
88  float dPhidRL1SeedP = initDRdPhi;
89  getL1MuonVariables(global_p, l1Muons, dR2dRL1SeedP, dPhidRL1SeedP);
90 
91  var[kDRdRL1SeedP] = std::sqrt(dR2dRL1SeedP);
92  var[kDPhidRL1SeedP] = dPhidRL1SeedP;
93 
94  if (!isFromL1_) {
95  float dR2dRL2SeedP = initDRdPhi;
96  float dPhidRL2SeedP = initDRdPhi;
97  getL2MuonVariables(global_p, l2Muons, dR2dRL2SeedP, dPhidRL2SeedP);
98 
99  var[kDRdRL2SeedP] = std::sqrt(dR2dRL2SeedP);
100  var[kDPhidRL2SeedP] = dPhidRL2SeedP;
101  }
102 
103  for (int iv = 0; iv < kLast; ++iv) {
104  var[iv] = (var[iv] - scale_mean_.at(iv)) / scale_std_.at(iv);
105  }
106 
107  return gbrForest_->GetResponse(var);
108 }
const bool isFromL1_
T sqrt(T t)
Definition: SSEVec.h:19
std::unique_ptr< const GBRForest > gbrForest_
const std::vector< double > scale_mean_
void getL2MuonVariables(const GlobalVector &, const reco::RecoChargedCandidateCollection &, float &, float &) const
void getL1MuonVariables(const GlobalVector &, const l1t::MuonBxCollection &, float &, float &) const
const std::vector< double > scale_std_

◆ getL1MuonVariables()

void SeedMvaEstimator::getL1MuonVariables ( const GlobalVector global_p,
const l1t::MuonBxCollection l1Muons,
float &  dR2dRL1SeedP,
float &  dPhidRL1SeedP 
) const
private

Definition at line 38 of file SeedMvaEstimator.cc.

References BXVector< T >::begin(), reco::deltaPhi(), reco::deltaR2(), BXVector< T >::end(), PV3DBase< T, PVType, FrameType >::eta(), BXVector< T >::getFirstBX(), BXVector< T >::getLastBX(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, minL1Qual_, and PV3DBase< T, PVType, FrameType >::phi().

Referenced by computeMva().

41  {
42  for (int ibx = l1Muons.getFirstBX(); ibx <= l1Muons.getLastBX(); ++ibx) {
43  if (ibx != 0)
44  continue; // -- only take when ibx == 0 -- //
45 
46  for (auto it = l1Muons.begin(ibx); it != l1Muons.end(ibx); it++) {
47  if (it->hwQual() < minL1Qual_)
48  continue;
49 
50  float dR2tmp = reco::deltaR2(it->etaAtVtx(), it->phiAtVtx(), global_p.eta(), global_p.phi());
51  if (dR2tmp < dR2dRL1SeedP) {
52  dR2dRL1SeedP = dR2tmp;
53  dPhidRL1SeedP = reco::deltaPhi(it->phiAtVtx(), global_p.phi());
54  }
55  }
56  }
57 }
int getLastBX() const
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
int getFirstBX() const
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
T eta() const
Definition: PV3DBase.h:73
const_iterator begin(int bx) const
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
const_iterator end(int bx) const

◆ getL2MuonVariables()

void SeedMvaEstimator::getL2MuonVariables ( const GlobalVector global_p,
const reco::RecoChargedCandidateCollection l2Muons,
float &  dR2dRL2SeedP,
float &  dPhidRL2SeedP 
) const
private

Definition at line 59 of file SeedMvaEstimator.cc.

References reco::deltaPhi(), reco::deltaR2(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, and PV3DBase< T, PVType, FrameType >::phi().

Referenced by computeMva().

62  {
63  for (auto it = l2Muons.begin(); it != l2Muons.end(); it++) {
64  float dR2tmp = reco::deltaR2(*it, global_p);
65  if (dR2tmp < dR2dRL2SeedP) {
66  dR2dRL2SeedP = dR2tmp;
67  dPhidRL2SeedP = reco::deltaPhi(it->phi(), global_p.phi());
68  }
69  }
70 }
constexpr double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:26
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16

Member Data Documentation

◆ gbrForest_

std::unique_ptr<const GBRForest> SeedMvaEstimator::gbrForest_
private

Definition at line 38 of file SeedMvaEstimator.h.

Referenced by computeMva(), and SeedMvaEstimator().

◆ isFromL1_

const bool SeedMvaEstimator::isFromL1_
private

Definition at line 41 of file SeedMvaEstimator.h.

Referenced by computeMva().

◆ minL1Qual_

const int SeedMvaEstimator::minL1Qual_
private

Definition at line 42 of file SeedMvaEstimator.h.

Referenced by getL1MuonVariables().

◆ scale_mean_

const std::vector<double> SeedMvaEstimator::scale_mean_
private

Definition at line 39 of file SeedMvaEstimator.h.

Referenced by computeMva().

◆ scale_std_

const std::vector<double> SeedMvaEstimator::scale_std_
private

Definition at line 40 of file SeedMvaEstimator.h.

Referenced by computeMva().