Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:11:15

0001 #include "FastSimulation/Event/interface/KineParticleFilter.h"
0002 #include "CommonTools/BaseParticlePropagator/interface/RawParticle.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 
0005 KineParticleFilter::KineParticleFilter(const edm::ParameterSet& cfg) {
0006   // Charged particles must have pt greater than chargedPtMin [GeV]
0007   double chargedPtMin = cfg.getParameter<double>("chargedPtMin");
0008   chargedPtMin2 = chargedPtMin * chargedPtMin;
0009 
0010   // Particles must have energy greater than EMin [GeV]
0011   EMin = cfg.getParameter<double>("EMin");
0012 
0013   // Allow *ALL* protons with energy > protonEMin
0014   protonEMin = cfg.getParameter<double>("protonEMin");
0015 
0016   // Particles must have abs(eta) < etaMax (if close enough to 0,0,0)
0017   double etaMax = cfg.getParameter<double>("etaMax");
0018   cos2ThetaMax = (std::exp(2. * etaMax) - 1.) / (std::exp(2. * etaMax) + 1.);
0019   cos2ThetaMax *= cos2ThetaMax;
0020 
0021   // Particles must have vertex inside the volume enclosed by ECAL
0022   double vertexRMax = cfg.getParameter<double>("rMax");
0023   vertexRMax2 = vertexRMax * vertexRMax;
0024   vertexZMax = cfg.getParameter<double>("zMax");
0025 }
0026 
0027 bool KineParticleFilter::acceptParticle(const RawParticle& particle) const {
0028   int pId = abs(particle.pid());
0029 
0030   // skipp invisible particles
0031   if (pId == 12 || pId == 14 || pId == 16 || pId == 1000022) {
0032     return false;
0033   }
0034 
0035   // keep all high-energy protons
0036   else if (pId == 2212 && particle.E() >= protonEMin) {
0037     return acceptVertex(particle.vertex());
0038   }
0039 
0040   // cut on the energy
0041   else if (particle.E() < EMin) {
0042     return false;
0043   }
0044 
0045   // cut on pt of charged particles
0046   else if (particle.charge() != 0 && particle.Perp2() < chargedPtMin2) {
0047     return false;
0048   }
0049 
0050   // cut on eta if the origin vertex is close to the beam
0051   else if (particle.vertex().Perp2() < 25. && particle.cos2Theta() > cos2ThetaMax) {
0052     return false;
0053   }
0054 
0055   // particles must have vertex in volume enclosed by ECAL
0056   return acceptVertex(particle.vertex());
0057 }
0058 
0059 bool KineParticleFilter::acceptVertex(const XYZTLorentzVector& vertex) const {
0060   return vertex.Perp2() < vertexRMax2 && fabs(vertex.Z()) < vertexZMax;
0061 }