CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
fastsim::ParticleFilter Class Reference

(Kinematic) cuts on the particles that are propagated. More...

#include <ParticleFilter.h>

Public Member Functions

bool accepts (const Particle &particle) const
 Check all if all criteria are fullfilled. More...
 
bool acceptsEn (const Particle &particle) const
 Kinematic cuts on the particle. More...
 
bool acceptsVtx (const math::XYZTLorentzVector &originVertexPosition) const
 Vertex within tracker volume. More...
 
 ParticleFilter (const edm::ParameterSet &cfg)
 Default Constructor. More...
 

Private Attributes

double chargedPtMin2_
 Minimum pT^2 of a charged particle. More...
 
double cos2ThetaMax_
 Particles must have abs(eta) < etaMax if close to beampipe. More...
 
double EMin_
 Minimum energy of a particle. More...
 
double protonEMin_
 Allow ALL protons with energy > protonEMin. More...
 
std::vector< int > skipParticles_
 List of invisible particles (neutrinos are excluded by default) More...
 
double vertexRMax2_
 Radius^2 of tracker volume. More...
 
double vertexZMax_
 Z of tracker volume. More...
 

Detailed Description

(Kinematic) cuts on the particles that are propagated.

All other particles are skipped.

Definition at line 30 of file ParticleFilter.h.

Constructor & Destructor Documentation

fastsim::ParticleFilter::ParticleFilter ( const edm::ParameterSet cfg)

Default Constructor.

Definition at line 6 of file ParticleFilter.cc.

References particleFlowSimParticle_cfi::chargedPtMin, chargedPtMin2_, cos2ThetaMax_, EMin_, ALCARECOTkAlBeamHalo_cff::etaMax, edm::ParameterSet::getParameter(), protonEMin_, skipParticles_, vertexRMax2_, and vertexZMax_.

7 {
8  // Charged particles must have pt greater than chargedPtMin [GeV]
9  double chargedPtMin = cfg.getParameter<double>("chargedPtMin");
10  chargedPtMin2_ = chargedPtMin*chargedPtMin;
11 
12  // (All) particles must have energy greater than EMin [GeV]
13  EMin_ = cfg.getParameter<double>("EMin");
14 
15  // Allow *ALL* protons with energy > protonEMin
16  protonEMin_ = cfg.getParameter<double>("protonEMin");
17 
18  // List of invisible particles if extension needed
19  // Predefined: Neutrinos, Neutralino_1
20  skipParticles_ = cfg.getParameter<std::vector<int>>("invisibleParticles");
21 
22  // Particles must have abs(eta) < etaMax (if close enough to 0,0,0)
23  double etaMax = cfg.getParameter<double>("etaMax");
24  cos2ThetaMax_ = (vdt::fast_exp(2.*etaMax)-1.) / (vdt::fast_exp(2.*etaMax)+1.);
25  cos2ThetaMax_ *= cos2ThetaMax_;
26 
27  // Particles must have vertex inside the tracker
28  vertexRMax2_ = 129.0*129.0;
29  vertexZMax_ = 303.353;
30 }
T getParameter(std::string const &) const
double cos2ThetaMax_
Particles must have abs(eta) < etaMax if close to beampipe.
double EMin_
Minimum energy of a particle.
std::vector< int > skipParticles_
List of invisible particles (neutrinos are excluded by default)
double protonEMin_
Allow ALL protons with energy > protonEMin.
double vertexRMax2_
Radius^2 of tracker volume.
double vertexZMax_
Z of tracker volume.
double chargedPtMin2_
Minimum pT^2 of a charged particle.

Member Function Documentation

bool fastsim::ParticleFilter::accepts ( const Particle particle) const

Check all if all criteria are fullfilled.

  • Particle is invisible (neutrinos by default, list can be extended)
  • Kinematic cuts (calls acceptsEN(...))
  • Vertex within tracker volume (calls acceptsVtx(...))
    See also
    acceptsEn(const Particle & particle)
    acceptsVtx(const math::XYZTLorentzVector & originVertexPosition)

Definition at line 32 of file ParticleFilter.cc.

References funct::abs(), acceptsEn(), acceptsVtx(), cos2ThetaMax_, fastsim::Particle::momentum(), fastsim::Particle::pdgId(), fastsim::Particle::position(), protonEMin_, and skipParticles_.

Referenced by fastsim::ParticleManager::nextParticle().

33 {
34  int absPdgId = abs(particle.pdgId());
35 
36  // skip invisible particles
37  if(absPdgId == 12 || absPdgId == 14 || absPdgId == 16 || absPdgId == 1000022)
38  {
39  return false;
40  }
41  // keep all high-energy protons
42  else if(absPdgId == 2212 && particle.momentum().E() >= protonEMin_)
43  {
44  return true;
45  }
46 
47  // cut on eta if the origin vertex is close to the beam
48  else if(particle.position().Perp2() < 25. && particle.momentum().Pz()*particle.momentum().Pz()/particle.momentum().P2() > cos2ThetaMax_)
49  {
50  return false;
51  }
52 
53  // possible to extend list of invisible particles
54  for(unsigned InvIdx = 0; InvIdx < skipParticles_.size(); InvIdx++){
55  if(absPdgId == abs(skipParticles_.at(InvIdx)))
56  {
57  return false;
58  }
59  }
60 
61  // particles must have vertex in volume of tracker
62  return acceptsVtx(particle.position()) && acceptsEn(particle);
63  //return (acceptsVtx(particle.position()) || particle.momentum().Pz()*particle.momentum().Pz()/particle.momentum().P2() > (vdt::fast_exp(2.*3.0)-1.) / (vdt::fast_exp(2.*3.0)+1.)*(vdt::fast_exp(2.*3.0)-1.) / (vdt::fast_exp(2.*3.0)+1.)) && acceptsEn(particle);
64 }
double cos2ThetaMax_
Particles must have abs(eta) < etaMax if close to beampipe.
std::vector< int > skipParticles_
List of invisible particles (neutrinos are excluded by default)
double protonEMin_
Allow ALL protons with energy > protonEMin.
bool acceptsEn(const Particle &particle) const
Kinematic cuts on the particle.
double Pz() const
Definition: Particle.h:61
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double E() const
Definition: Particle.h:63
bool acceptsVtx(const math::XYZTLorentzVector &originVertexPosition) const
Vertex within tracker volume.
bool fastsim::ParticleFilter::acceptsEn ( const Particle particle) const

Kinematic cuts on the particle.

Definition at line 66 of file ParticleFilter.cc.

References funct::abs(), fastsim::Particle::charge(), chargedPtMin2_, EMin_, fastsim::Particle::momentum(), fastsim::Particle::pdgId(), and protonEMin_.

Referenced by accepts().

67 {
68  int absPdgId = abs(particle.pdgId());
69 
70  // keep all high-energy protons
71  if(absPdgId == 2212 && particle.momentum().E() >= protonEMin_)
72  {
73  return true;
74  }
75 
76  // cut on the energy
77  else if(particle.momentum().E() < EMin_)
78  {
79  return false;
80  }
81 
82  // cut on pt of charged particles
83  else if(particle.charge()!=0 && particle.momentum().Perp2()<chargedPtMin2_)
84  {
85  return false;
86  }
87 
88  return true;
89 }
double EMin_
Minimum energy of a particle.
double protonEMin_
Allow ALL protons with energy > protonEMin.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double chargedPtMin2_
Minimum pT^2 of a charged particle.
double E() const
Definition: Particle.h:63
bool fastsim::ParticleFilter::acceptsVtx ( const math::XYZTLorentzVector originVertexPosition) const

Vertex within tracker volume.

Parameters
originVertexPositionPosition of origin vertex.

Definition at line 92 of file ParticleFilter.cc.

References vertexRMax2_, and vertexZMax_.

Referenced by accepts(), and fastsim::ParticleManager::addSecondaries().

93 {
94  // origin vertex is within the tracker volume
95  return (originVertex.Perp2() < vertexRMax2_ && fabs(originVertex.Z()) < vertexZMax_);
96 }
double vertexRMax2_
Radius^2 of tracker volume.
double vertexZMax_
Z of tracker volume.

Member Data Documentation

double fastsim::ParticleFilter::chargedPtMin2_
private

Minimum pT^2 of a charged particle.

Definition at line 56 of file ParticleFilter.h.

Referenced by acceptsEn(), and ParticleFilter().

double fastsim::ParticleFilter::cos2ThetaMax_
private

Particles must have abs(eta) < etaMax if close to beampipe.

Definition at line 59 of file ParticleFilter.h.

Referenced by accepts(), and ParticleFilter().

double fastsim::ParticleFilter::EMin_
private

Minimum energy of a particle.

Definition at line 57 of file ParticleFilter.h.

Referenced by acceptsEn(), and ParticleFilter().

double fastsim::ParticleFilter::protonEMin_
private

Allow ALL protons with energy > protonEMin.

Definition at line 58 of file ParticleFilter.h.

Referenced by accepts(), acceptsEn(), and ParticleFilter().

std::vector<int> fastsim::ParticleFilter::skipParticles_
private

List of invisible particles (neutrinos are excluded by default)

Definition at line 62 of file ParticleFilter.h.

Referenced by accepts(), and ParticleFilter().

double fastsim::ParticleFilter::vertexRMax2_
private

Radius^2 of tracker volume.

Definition at line 60 of file ParticleFilter.h.

Referenced by acceptsVtx(), and ParticleFilter().

double fastsim::ParticleFilter::vertexZMax_
private

Z of tracker volume.

Definition at line 61 of file ParticleFilter.h.

Referenced by acceptsVtx(), and ParticleFilter().