CMS 3D CMS Logo

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

#include <GenParticleCustomSelector.h>

Public Member Functions

 GenParticleCustomSelector ()
 
 GenParticleCustomSelector (double ptMin, double minRapidity, double maxRapidity, double tip, double lip, bool chargedOnly, int status, const std::vector< int > &pdgId=std::vector< int >(), bool invertRapidityCut=false, double minPhi=-3.2, double maxPhi=3.2)
 
bool operator() (const reco::GenParticle &tp) const
 Operator() performs the selection: e.g. if (tPSelector(tp)) {...}. More...
 

Private Attributes

bool chargedOnly_
 
bool invertRapidityCut_
 
double lip_
 
double maxRapidity_
 
float meanPhi_
 
double minRapidity_
 
std::vector< int > pdgId_
 
double ptMin_
 
float rangePhi_
 
int status_
 
double tip_
 

Detailed Description

Definition at line 11 of file GenParticleCustomSelector.h.

Constructor & Destructor Documentation

◆ GenParticleCustomSelector() [1/2]

GenParticleCustomSelector::GenParticleCustomSelector ( )
inline

Definition at line 13 of file GenParticleCustomSelector.h.

13 {}

◆ GenParticleCustomSelector() [2/2]

GenParticleCustomSelector::GenParticleCustomSelector ( double  ptMin,
double  minRapidity,
double  maxRapidity,
double  tip,
double  lip,
bool  chargedOnly,
int  status,
const std::vector< int > &  pdgId = std::vector<int>(),
bool  invertRapidityCut = false,
double  minPhi = -3.2,
double  maxPhi = 3.2 
)
inline

Definition at line 14 of file GenParticleCustomSelector.h.

References Exception, M_PI, HLT_2024v12_cff::maxPhi, and HLT_2024v12_cff::minPhi.

25  : ptMin_(ptMin),
28  meanPhi_((minPhi + maxPhi) / 2.),
29  rangePhi_((maxPhi - minPhi) / 2.),
30  tip_(tip),
31  lip_(lip),
33  status_(status),
34  pdgId_(pdgId),
36  if (minPhi >= maxPhi) {
37  throw cms::Exception("Configuration")
38  << "GenParticleCustomSelector: minPhi (" << minPhi << ") must be smaller than maxPhi (" << maxPhi
39  << "). The range is constructed from minPhi to maxPhi around their "
40  "average.";
41  }
42  if (minPhi >= M_PI) {
43  throw cms::Exception("Configuration") << "GenParticleCustomSelector: minPhi (" << minPhi
44  << ") must be smaller than PI. The range is constructed from minPhi "
45  "to maxPhi around their average.";
46  }
47  if (maxPhi <= -M_PI) {
48  throw cms::Exception("Configuration") << "GenParticleCustomSelector: maxPhi (" << maxPhi
49  << ") must be larger than -PI. The range is constructed from minPhi "
50  "to maxPhi around their average.";
51  }
52  }
constexpr float ptMin
#define M_PI

Member Function Documentation

◆ operator()()

bool GenParticleCustomSelector::operator() ( const reco::GenParticle tp) const
inline

Operator() performs the selection: e.g. if (tPSelector(tp)) {...}.

Definition at line 55 of file GenParticleCustomSelector.h.

References chargedOnly_, SiPixelRawToDigiRegional_cfi::deltaPhi, PVValHelper::eta, invertRapidityCut_, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, lip_, maxRapidity_, meanPhi_, minRapidity_, AlCaHLTBitMon_ParallelJobs::p, pdgId_, DiDispStaMuonMonitor_cfi::pt, ptMin_, rangePhi_, mathSSE::sqrt(), status_, tip_, and cmsswSequenceInfo::tp.

55  {
56  if (chargedOnly_ && tp.charge() == 0)
57  return false; //select only if charge!=0
58  bool testId = false;
59  unsigned int idSize = pdgId_.size();
60  if (idSize == 0)
61  testId = true;
62  else
63  for (unsigned int it = 0; it != idSize; ++it) {
64  if (tp.pdgId() == pdgId_[it])
65  testId = true;
66  }
67 
68  auto etaOk = [&](const reco::GenParticle& p) -> bool {
69  float eta = p.eta();
70  if (!invertRapidityCut_)
71  return (eta >= minRapidity_) && (eta <= maxRapidity_);
72  else
73  return (eta < minRapidity_ || eta > maxRapidity_);
74  };
75  auto phiOk = [&](const reco::GenParticle& p) {
76  float dphi = deltaPhi(atan2f(p.py(), p.px()), meanPhi_);
77  return dphi >= -rangePhi_ && dphi <= rangePhi_;
78  };
79  auto ptOk = [&](const reco::GenParticle& p) {
80  double pt = p.pt();
81  return pt >= ptMin_;
82  };
83 
84  return (ptOk(tp) && etaOk(tp) && phiOk(tp) && sqrt(tp.vertex().perp2()) <= tip_ && fabs(tp.vertex().z()) <= lip_ &&
85  tp.status() == status_ && testId);
86  }
T sqrt(T t)
Definition: SSEVec.h:19

Member Data Documentation

◆ chargedOnly_

bool GenParticleCustomSelector::chargedOnly_
private

Definition at line 96 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ invertRapidityCut_

bool GenParticleCustomSelector::invertRapidityCut_
private

Definition at line 99 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ lip_

double GenParticleCustomSelector::lip_
private

Definition at line 95 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ maxRapidity_

double GenParticleCustomSelector::maxRapidity_
private

Definition at line 91 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ meanPhi_

float GenParticleCustomSelector::meanPhi_
private

Definition at line 92 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ minRapidity_

double GenParticleCustomSelector::minRapidity_
private

Definition at line 90 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ pdgId_

std::vector<int> GenParticleCustomSelector::pdgId_
private

Definition at line 98 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ ptMin_

double GenParticleCustomSelector::ptMin_
private

Definition at line 89 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ rangePhi_

float GenParticleCustomSelector::rangePhi_
private

Definition at line 93 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ status_

int GenParticleCustomSelector::status_
private

Definition at line 97 of file GenParticleCustomSelector.h.

Referenced by operator()().

◆ tip_

double GenParticleCustomSelector::tip_
private

Definition at line 94 of file GenParticleCustomSelector.h.

Referenced by operator()().