hej is hosted by Hepforge, IPPP Durham
HEJ 2 2.0
High energy resummation for hadron colliders
Loading...
Searching...
No Matches
Particle.hh
Go to the documentation of this file.
1
10#pragma once
11
12#include "fastjet/PseudoJet.hh"
13
14#include "HEJ/PDG_codes.hh"
15
16namespace HEJ {
17
19 struct Particle {
23 fastjet::PseudoJet p;
24
26 double rapidity() const{
27 return p.rapidity();
28 }
30 double perp() const{
31 return p.perp();
32 }
34 double px() const{
35 return p.px();
36 }
38 double py() const{
39 return p.py();
40 }
42 double pz() const{
43 return p.pz();
44 }
46 double E() const{
47 return p.E();
48 }
50 double m() const{
51 return p.m();
52 }
53 };
54
56
63 template<class FourVector>
64 bool operator()(FourVector const & p1, FourVector const & p2){
65 return p1.rapidity() < p2.rapidity();
66 }
67 };
68
70
76 struct pz_less{
77 template<class FourVector>
78 bool operator()(FourVector const & p1, FourVector const & p2){
79 return p1.pz() < p2.pz();
80 }
81 };
82
83
85 inline
86 std::vector<fastjet::PseudoJet> to_PseudoJet(
87 std::vector<Particle> const & v
88 ){
89 std::vector<fastjet::PseudoJet> result;
90 for(auto && sp: v) result.emplace_back(sp.p);
91 return result;
92 }
93
95 inline
96 bool is_parton(Particle const & p){
97 return is_parton(p.type);
98 }
99
101 inline bool is_AWZH_boson(Particle const & particle){
102 return is_AWZH_boson(particle.type);
103 }
104
106 inline
107 std::vector<Particle> filter_partons(
108 std::vector<Particle> const & v
109 ){
110 std::vector<Particle> result;
111 result.reserve(v.size());
112 std::copy_if(
113 begin(v), end(v), std::back_inserter(result),
114 [](Particle const & p){ return is_parton(p); }
115 );
116 return result;
117 }
118}
Contains the Particle IDs of all relevant SM particles.
ParticleID
The possible particle identities. We use PDG IDs as standard.
Definition: PDG_codes.hh:23
Main HEJ 2 Namespace.
Definition: mainpage.dox:1
std::vector< fastjet::PseudoJet > to_PseudoJet(std::vector< Particle > const &v)
Convert a vector of Particles to a vector of particle momenta.
Definition: Particle.hh:86
std::vector< Particle > filter_partons(std::vector< Particle > const &v)
Extract all partons from a vector of particles.
Definition: Particle.hh:107
bool is_parton(Particle const &p)
Check if a particle is a parton, i.e. quark, antiquark, or gluon.
Definition: Particle.hh:96
bool is_AWZH_boson(Particle const &particle)
Check if a particle is a photon, W, Z, or Higgs boson.
Definition: Particle.hh:101
Class representing a particle.
Definition: Particle.hh:19
double rapidity() const
get rapidity
Definition: Particle.hh:26
double py() const
get momentum in y direction
Definition: Particle.hh:38
fastjet::PseudoJet p
particle momentum
Definition: Particle.hh:23
double px() const
get momentum in x direction
Definition: Particle.hh:34
double E() const
get energy
Definition: Particle.hh:46
ParticleID type
particle type
Definition: Particle.hh:21
double pz() const
get momentum in z direction
Definition: Particle.hh:42
double m() const
get mass
Definition: Particle.hh:50
double perp() const
get transverse momentum
Definition: Particle.hh:30
Functor to compare momenta in z direction.
Definition: Particle.hh:76
bool operator()(FourVector const &p1, FourVector const &p2)
Definition: Particle.hh:78
Functor to compare rapidities.
Definition: Particle.hh:62
bool operator()(FourVector const &p1, FourVector const &p2)
Definition: Particle.hh:64