CMS 3D CMS Logo

KineParticleFilter.cc
Go to the documentation of this file.
4 
6 {
7  // Charged particles must have pt greater than chargedPtMin [GeV]
8  double chargedPtMin = cfg.getParameter<double>("chargedPtMin");
9  chargedPtMin2 = chargedPtMin*chargedPtMin;
10 
11  // Particles must have energy greater than EMin [GeV]
12  EMin = cfg.getParameter<double>("EMin");
13 
14  // Allow *ALL* protons with energy > protonEMin
15  protonEMin = cfg.getParameter<double>("protonEMin");
16 
17  // Particles must have abs(eta) < etaMax (if close enough to 0,0,0)
18  double etaMax = cfg.getParameter<double>("etaMax");
19  cos2ThetaMax = (std::exp(2.*etaMax)-1.) / (std::exp(2.*etaMax)+1.);
20  cos2ThetaMax *= cos2ThetaMax;
21 
22  // Particles must have vertex inside the volume enclosed by ECAL
23  double vertexRMax = cfg.getParameter<double>("rMax");
24  vertexRMax2 = vertexRMax*vertexRMax;
25  vertexZMax = cfg.getParameter<double>("zMax");
26 
27 }
28 
29 bool KineParticleFilter::acceptParticle(const RawParticle & particle) const
30 {
31 
32  int pId = abs(particle.pid());
33 
34  // skipp invisible particles
35  if(pId == 12 || pId == 14 || pId == 16 || pId == 1000022)
36  {
37  return false;
38  }
39 
40  // keep all high-energy protons
41  else if(pId == 2212 && particle.E() >= protonEMin)
42  {
43  return true;
44  }
45 
46  // cut on the energy
47  else if( particle.E() < EMin)
48  {
49  return false;
50  }
51 
52  // cut on pt of charged particles
53  else if( particle.charge()!=0 && particle.Perp2()<chargedPtMin2)
54  {
55  return false;
56  }
57 
58  // cut on eta if the origin vertex is close to the beam
59  else if( particle.vertex().Perp2() < 25. && particle.cos2Theta() > cos2ThetaMax)
60  {
61  return false;
62  }
63 
64 
65 
66  // particles must have vertex in volume enclosed by ECAL
67  return acceptVertex(particle.vertex());
68 }
69 
70 
72 {
73  return vertex.Perp2() < vertexRMax2 && fabs(vertex.Z()) < vertexZMax;
74 }
T getParameter(std::string const &) const
bool acceptParticle(const RawParticle &p) const
double Perp2() const
perpendicular momentum squared
Definition: RawParticle.h:330
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:296
bool acceptVertex(const math::XYZTLorentzVector &p) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double charge() const
get the MEASURED charge
Definition: RawParticle.h:313
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:339
KineParticleFilter(const edm::ParameterSet &kine)
double E() const
energy of the momentum
Definition: RawParticle.h:325
double cos2Theta() const
Cos**2(theta) is faster to determine than eta.
Definition: RawParticle.h:299
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:27