CMS 3D CMS Logo

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

#include <CaloParticleSelector.h>

Public Member Functions

 CaloParticleSelector ()
 
 CaloParticleSelector (double ptMin, double ptMax, double minRapidity, double maxRapidity, double lip, double tip, int minHit, unsigned int maxSimClusters, bool signalOnly, bool intimeOnly, bool chargedOnly, bool stableOnly, bool notConvertedOnly, const std::vector< int > &pdgId=std::vector< int >(), double minPhi=-3.2, double maxPhi=3.2)
 
bool operator() (const CaloParticle &cp, std::vector< SimVertex > const &simVertices) const
 

Private Attributes

bool chargedOnly_
 
bool intimeOnly_
 
double lip_
 
float maxRapidity_
 
unsigned int maxSimClusters_
 
float meanPhi_
 
int minHit_
 
float minRapidity_
 
bool notConvertedOnly_
 
std::vector< int > pdgId_
 
double ptMax2_
 
double ptMin2_
 
float rangePhi_
 
bool signalOnly_
 
bool stableOnly_
 
double tip2_
 

Detailed Description

Definition at line 10 of file CaloParticleSelector.h.

Constructor & Destructor Documentation

◆ CaloParticleSelector() [1/2]

CaloParticleSelector::CaloParticleSelector ( )
inline

Definition at line 12 of file CaloParticleSelector.h.

12 {}

◆ CaloParticleSelector() [2/2]

CaloParticleSelector::CaloParticleSelector ( double  ptMin,
double  ptMax,
double  minRapidity,
double  maxRapidity,
double  lip,
double  tip,
int  minHit,
unsigned int  maxSimClusters,
bool  signalOnly,
bool  intimeOnly,
bool  chargedOnly,
bool  stableOnly,
bool  notConvertedOnly,
const std::vector< int > &  pdgId = std::vector<int>(),
double  minPhi = -3.2,
double  maxPhi = 3.2 
)
inline

Definition at line 13 of file CaloParticleSelector.h.

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

29  : ptMin2_(ptMin * ptMin),
30  ptMax2_(ptMax * ptMax),
33  lip_(lip),
34  tip2_(tip * tip),
35  meanPhi_((minPhi + maxPhi) / 2.),
36  rangePhi_((maxPhi - minPhi) / 2.),
37  minHit_(minHit),
38  maxSimClusters_(maxSimClusters),
43  notConvertedOnly_(notConvertedOnly),
44  pdgId_(pdgId) {
45  if (minPhi >= maxPhi) {
46  throw cms::Exception("Configuration")
47  << "CaloParticleSelector: minPhi (" << minPhi << ") must be smaller than maxPhi (" << maxPhi
48  << "). The range is constructed from minPhi to maxPhi around their average.";
49  }
50  if (minPhi >= M_PI) {
51  throw cms::Exception("Configuration")
52  << "CaloParticleSelector: minPhi (" << minPhi
53  << ") must be smaller than PI. The range is constructed from minPhi to maxPhi around their average.";
54  }
55  if (maxPhi <= -M_PI) {
56  throw cms::Exception("Configuration")
57  << "CaloParticleSelector: maxPhi (" << maxPhi
58  << ") must be larger than -PI. The range is constructed from minPhi to maxPhi around their average.";
59  }
60  }
constexpr float ptMin
#define M_PI
std::vector< int > pdgId_

Member Function Documentation

◆ operator()()

bool CaloParticleSelector::operator() ( const CaloParticle cp,
std::vector< SimVertex > const &  simVertices 
) const
inline

Definition at line 64 of file CaloParticleSelector.h.

References funct::abs(), chargedOnly_, SiPixelRawToDigiRegional_cfi::deltaPhi, PVValHelper::eta, intimeOnly_, dqmiolumiharvest::j, maxRapidity_, maxSimClusters_, meanPhi_, minRapidity_, notConvertedOnly_, AlCaHLTBitMon_ParallelJobs::p, EgammaValidation_cff::pdgid, pdgId_, ntuple::pdgids, HLT_2024v11_cff::pt2, ptMax2_, ptMin2_, rangePhi_, signalOnly_, and stableOnly_.

64  {
65  // signal only means no PU particles
66  if (signalOnly_ && !(cp.eventId().bunchCrossing() == 0 && cp.eventId().event() == 0))
67  return false;
68  // intime only means no OOT PU particles
69  if (intimeOnly_ && !(cp.eventId().bunchCrossing() == 0))
70  return false;
71 
72  if (cp.simClusters().size() > maxSimClusters_)
73  return false;
74 
75  auto pdgid = cp.pdgId();
76  if (!pdgId_.empty()) {
77  bool testId = false;
78  for (auto id : pdgId_) {
79  if (id == pdgid) {
80  testId = true;
81  break;
82  }
83  }
84  if (!testId)
85  return false;
86  }
87 
88  if (chargedOnly_ && cp.charge() == 0)
89  return false; //select only if charge!=0
90 
91  // select only stable particles
92  if (stableOnly_) {
93  for (CaloParticle::genp_iterator j = cp.genParticle_begin(); j != cp.genParticle_end(); ++j) {
94  if (j->get() == nullptr || j->get()->status() != 1) {
95  return false;
96  }
97  }
98 
99  // test for remaining unstabled due to lack of genparticle pointer
100  std::vector<int> pdgids{11, 13, 211, 321, 2212, 3112, 3222, 3312, 3334};
101  if (cp.status() == -99 && (!std::binary_search(pdgids.begin(), pdgids.end(), std::abs(pdgid)))) {
102  return false;
103  }
104  }
105 
106  // select only particles which did not convert/decay before the calorimeter
107  // in case of electrons, bremsstrahlung is usually the cause, thus this selection is skipped
108  if (std::abs(pdgid) != 11) {
109  if (notConvertedOnly_) {
110  if (cp.g4Tracks()[0].getPositionAtBoundary() == math::XYZTLorentzVectorF(0, 0, 0, 0)) {
111  return false;
112  }
113  if (cp.g4Tracks()[0].getMomentumAtBoundary() == math::XYZTLorentzVectorF(0, 0, 0, 0)) {
114  return false;
115  }
116  }
117  }
118 
119  auto etaOk = [&](const CaloParticle& p) -> bool {
120  float eta = etaFromXYZ(p.px(), p.py(), p.pz());
121  return (eta >= minRapidity_) & (eta <= maxRapidity_);
122  };
123  auto phiOk = [&](const CaloParticle& p) {
124  float dphi = deltaPhi(atan2f(p.py(), p.px()), meanPhi_);
125  return dphi >= -rangePhi_ && dphi <= rangePhi_;
126  };
127  auto ptOk = [&](const CaloParticle& p) {
128  double pt2 = cp.p4().perp2();
129  return pt2 >= ptMin2_ && pt2 <= ptMax2_;
130  };
131 
132  return (ptOk(cp) && etaOk(cp) && phiOk(cp));
133  }
pdgids
Definition: ntuple.py:97
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > pdgId_
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< float > > XYZTLorentzVectorF
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:22

Member Data Documentation

◆ chargedOnly_

bool CaloParticleSelector::chargedOnly_
private

Definition at line 148 of file CaloParticleSelector.h.

Referenced by operator()().

◆ intimeOnly_

bool CaloParticleSelector::intimeOnly_
private

Definition at line 147 of file CaloParticleSelector.h.

Referenced by operator()().

◆ lip_

double CaloParticleSelector::lip_
private

Definition at line 140 of file CaloParticleSelector.h.

◆ maxRapidity_

float CaloParticleSelector::maxRapidity_
private

Definition at line 139 of file CaloParticleSelector.h.

Referenced by operator()().

◆ maxSimClusters_

unsigned int CaloParticleSelector::maxSimClusters_
private

Definition at line 145 of file CaloParticleSelector.h.

Referenced by operator()().

◆ meanPhi_

float CaloParticleSelector::meanPhi_
private

Definition at line 142 of file CaloParticleSelector.h.

Referenced by operator()().

◆ minHit_

int CaloParticleSelector::minHit_
private

Definition at line 144 of file CaloParticleSelector.h.

◆ minRapidity_

float CaloParticleSelector::minRapidity_
private

Definition at line 138 of file CaloParticleSelector.h.

Referenced by operator()().

◆ notConvertedOnly_

bool CaloParticleSelector::notConvertedOnly_
private

Definition at line 150 of file CaloParticleSelector.h.

Referenced by operator()().

◆ pdgId_

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

Definition at line 151 of file CaloParticleSelector.h.

Referenced by operator()().

◆ ptMax2_

double CaloParticleSelector::ptMax2_
private

Definition at line 137 of file CaloParticleSelector.h.

Referenced by operator()().

◆ ptMin2_

double CaloParticleSelector::ptMin2_
private

Definition at line 136 of file CaloParticleSelector.h.

Referenced by operator()().

◆ rangePhi_

float CaloParticleSelector::rangePhi_
private

Definition at line 143 of file CaloParticleSelector.h.

Referenced by operator()().

◆ signalOnly_

bool CaloParticleSelector::signalOnly_
private

Definition at line 146 of file CaloParticleSelector.h.

Referenced by operator()().

◆ stableOnly_

bool CaloParticleSelector::stableOnly_
private

Definition at line 149 of file CaloParticleSelector.h.

Referenced by operator()().

◆ tip2_

double CaloParticleSelector::tip2_
private

Definition at line 141 of file CaloParticleSelector.h.