FCCAnalyses
Loading...
Searching...
No Matches
MCParticle.h
Go to the documentation of this file.
1
2#ifndef MCPARTICLE_ANALYZERS_H
3#define MCPARTICLE_ANALYZERS_H
4
5#include <cmath>
6#include <vector>
7
8#include "ROOT/RVec.hxx"
9#include "TLorentzVector.h"
10#include "edm4hep/MCParticleData.h"
11#include "edm4hep/ParticleIDData.h"
12#include "edm4hep/Vector3d.h"
13#include "edm4hep/Vector3f.h"
14
15namespace FCCAnalyses{
16
23namespace MCParticle{
24
26 struct filter_pdgID {
27 filter_pdgID(int arg_pdgid, bool arg_abs);
28 int m_pdgid; //> Generator pdgid
29 bool m_abs;//> Use absolute value for pdgig
30 bool operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
31 };
32
34 struct sel_pt {
35 sel_pt(float arg_min_pt);
36 float m_min_pt = 20; //> transverse momentum threshold [GeV]
37 ROOT::VecOps::RVec<edm4hep::MCParticleData> operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
38 };
39
42 sel_genStatus(int arg_status);
43 int m_status = 1; //> Generator status
44 ROOT::VecOps::RVec<edm4hep::MCParticleData> operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
45 };
46
48 struct sel_pdgID {
49 sel_pdgID(int arg_pdg, bool arg_chargeconjugate);
50 int m_pdg = 13;
51 bool m_chargeconjugate = true;
52 ROOT::VecOps::RVec<edm4hep::MCParticleData> operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
53 };
54
56 struct get_tree{
57 get_tree(int arg_index);
58 int m_index; //> MC Particle index to build the tree from
59 ROOT::VecOps::RVec<int> operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in, ROOT::VecOps::RVec<int> ind);
60 };
61
63 struct get_decay {
64 get_decay(int arg_mother, int arg_daughters, bool arg_inf);
65 int m_mother = 0; //> mother pdg id
66 int m_daughters = 0;//> daughters pdg id
67 bool m_inf = false;//> boolean to check if the pdgid is below a value rather than equal
68 bool operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in, ROOT::VecOps::RVec<int> ind);
69 };
70
73 get_EventPrimaryVertex( int arg_genstatus );
74 int m_genstatus = 21; // Pythia8 code of the incoming particles of the hardest subprocess
75 TVector3 operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
76 };
77
81 int m_genstatus = 21; // Pythia8 code of the incoming particles of the hardest subprocess
82 TLorentzVector operator() (ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
83 };
84
85
88 get_indices( int pdg_mother, std::vector<int> pdg_daughters, bool stableDaughters, bool chargeConjugateMother, bool chargeConjugateDaughters, bool inclusiveDecay ) ;
90 std::vector<int> m_pdg_daughters;
95 ROOT::VecOps::RVec<int> operator() ( ROOT::VecOps::RVec<edm4hep::MCParticleData> in , ROOT::VecOps::RVec<int> ind);
96 };
97
100 get_indices_ExclusiveDecay( int pdg_mother, std::vector<int> pdg_daughters, bool stableDaughters, bool chargeConjugate ) ;
101 };
102
104 ROOT::VecOps::RVec<int> get_indices_MotherByIndex( int imother,
105 std::vector<int> m_pdg_daughters,
106 bool m_stableDaughters,
107 bool m_chargeConjugateDaughters,
108 bool m_inclusiveDecay,
109 ROOT::VecOps::RVec<edm4hep::MCParticleData> in ,
110 ROOT::VecOps::RVec<int> ind);
111
113 ROOT::VecOps::RVec<int> get_indices_ExclusiveDecay_MotherByIndex( int imother,
114 std::vector<int> m_pdg_daughters,
115 bool m_stableDaughters,
116 ROOT::VecOps::RVec<edm4hep::MCParticleData> in ,
117 ROOT::VecOps::RVec<int> ind);
118
119
121 ROOT::VecOps::RVec<int> get_parentid(ROOT::VecOps::RVec<int> mcind, ROOT::VecOps::RVec<edm4hep::MCParticleData> mc, ROOT::VecOps::RVec<int> parents);
122
124 ROOT::VecOps::RVec<float> get_time(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
125
127 ROOT::VecOps::RVec<float> get_pdg(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
128
130 ROOT::VecOps::RVec<float> get_genStatus(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
131
133 ROOT::VecOps::RVec<float> get_simStatus(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
134
136 ROOT::VecOps::RVec<edm4hep::Vector3d> get_vertex(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
137
139 ROOT::VecOps::RVec<float> get_vertex_x(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
140
142 ROOT::VecOps::RVec<float> get_vertex_y(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
143
145 ROOT::VecOps::RVec<float> get_vertex_z(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
146
148 ROOT::VecOps::RVec<edm4hep::Vector3d> get_endPoint(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
149
151 ROOT::VecOps::RVec<edm4hep::Vector3d> get_endPoint(ROOT::VecOps::RVec<edm4hep::MCParticleData> in, ROOT::VecOps::RVec<int> ind );
152
154 ROOT::VecOps::RVec<float> get_endPoint_x(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
155
157 ROOT::VecOps::RVec<float> get_endPoint_y(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
158
160 ROOT::VecOps::RVec<float> get_endPoint_z(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
161
163 ROOT::VecOps::RVec<float> get_pt(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
164
166 ROOT::VecOps::RVec<float> get_p(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
167
169 ROOT::VecOps::RVec<float> get_px(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
170
172 ROOT::VecOps::RVec<float> get_py(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
173
175 ROOT::VecOps::RVec<float> get_pz(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
176
178 ROOT::VecOps::RVec<float> get_eta(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
179
181 ROOT::VecOps::RVec<float> get_y(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
182
184 ROOT::VecOps::RVec<float> get_theta(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
185
187 ROOT::VecOps::RVec<float> get_phi(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
188
190 ROOT::VecOps::RVec<float> get_e(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
191
193 ROOT::VecOps::RVec<float> get_mass(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
194
196 ROOT::VecOps::RVec<float> get_charge(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
197
199 ROOT::VecOps::RVec<TLorentzVector> get_tlv(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
200
202 ROOT::VecOps::RVec<edm4hep::MCParticleData> mergeParticles(ROOT::VecOps::RVec<edm4hep::MCParticleData> x, ROOT::VecOps::RVec<edm4hep::MCParticleData> y);
203
205 int get_n(ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
206
208 ROOT::VecOps::RVec<float> AngleBetweenTwoMCParticles( ROOT::VecOps::RVec<edm4hep::MCParticleData> p1, ROOT::VecOps::RVec<edm4hep::MCParticleData> p2 );
209
211 std::vector<int> get_list_of_stable_particles_from_decay( int i, ROOT::VecOps::RVec<edm4hep::MCParticleData> in, ROOT::VecOps::RVec<int> ind) ;
212
214 std::vector<int> get_list_of_particles_from_decay( int i,
215 ROOT::VecOps::RVec<edm4hep::MCParticleData> in,
216 ROOT::VecOps::RVec<int> ind);
217
219 edm4hep::MCParticleData sel_byIndex( int idx, ROOT::VecOps::RVec<edm4hep::MCParticleData> in);
220
222 std::vector<int> list_of_stable_particles_from_decay( int i, ROOT::VecOps::RVec<edm4hep::MCParticleData> in, ROOT::VecOps::RVec<int> ind) ;
224 std::vector<int> list_of_particles_from_decay( int i,
225 ROOT::VecOps::RVec<edm4hep::MCParticleData> in,
226 ROOT::VecOps::RVec<int> ind);
227
228
230 int get_lepton_origin(int idx,
231 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &in,
232 const ROOT::VecOps::RVec<int> &ind);
233
234 int get_lepton_origin(const edm4hep::MCParticleData &p,
235 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &in,
236 const ROOT::VecOps::RVec<int> &ind);
237
238 ROOT::VecOps::RVec<int> get_leptons_origin(const ROOT::VecOps::RVec<edm4hep::MCParticleData> &particles,
239 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &in,
240 const ROOT::VecOps::RVec<int> &ind);
241
242
243}//end NS MCParticle
244
245}//end NS FCCAnalyses
246#endif
ROOT::VecOps::RVec< int > get_leptons_origin(const ROOT::VecOps::RVec< edm4hep::MCParticleData > &particles, const ROOT::VecOps::RVec< edm4hep::MCParticleData > &in, const ROOT::VecOps::RVec< int > &ind)
Definition MCParticle.cc:764
edm4hep::MCParticleData sel_byIndex(int idx, ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
returns one MCParticle selected by its index in the particle block
Definition MCParticle.cc:470
ROOT::VecOps::RVec< float > get_px(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the momenta of the input MCParticles
Definition MCParticle.cc:360
ROOT::VecOps::RVec< float > get_simStatus(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the simulation status of the input MCParticles
Definition MCParticle.cc:218
ROOT::VecOps::RVec< float > get_py(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the momenta of the input MCParticles
Definition MCParticle.cc:368
ROOT::VecOps::RVec< edm4hep::Vector3d > get_vertex(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the production vertex of the input MCParticles
Definition MCParticle.cc:226
ROOT::VecOps::RVec< float > get_phi(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the phi of the input MCParticles
Definition MCParticle.cc:330
int get_lepton_origin(int idx, const ROOT::VecOps::RVec< edm4hep::MCParticleData > &in, const ROOT::VecOps::RVec< int > &ind)
return the pdg ID of the parent of a lepton (pre-FSR)
Definition MCParticle.cc:755
ROOT::VecOps::RVec< float > get_eta(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the pseudo-rapidity of the input MCParticles
Definition MCParticle.cc:320
ROOT::VecOps::RVec< int > get_indices_MotherByIndex(int imother, std::vector< int > m_pdg_daughters, bool m_stableDaughters, bool m_chargeConjugateDaughters, bool m_inclusiveDecay, ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
return a list of indices that correspond to a given MC decay
Definition MCParticle.cc:565
ROOT::VecOps::RVec< float > get_theta(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the theta of the input MCParticles
Definition MCParticle.cc:402
ROOT::VecOps::RVec< float > get_pz(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the momenta of the input MCParticles
Definition MCParticle.cc:376
ROOT::VecOps::RVec< float > get_genStatus(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the generator status of the input MCParticles
Definition MCParticle.cc:210
std::vector< int > list_of_particles_from_decay(int i, ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
obsolete: should use get_list_of_particles_from_decay instead
Definition MCParticle.cc:553
ROOT::VecOps::RVec< TLorentzVector > get_tlv(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the TlorentzVector of the input MCParticles
Definition MCParticle.cc:412
ROOT::VecOps::RVec< float > get_endPoint_y(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the end point y of the input MCParticles
Definition MCParticle.cc:296
ROOT::VecOps::RVec< float > get_endPoint_z(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the z of the input MCParticles
Definition MCParticle.cc:304
ROOT::VecOps::RVec< int > get_indices_ExclusiveDecay_MotherByIndex(int imother, std::vector< int > m_pdg_daughters, bool m_stableDaughters, ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
a shorthand for get_indices_MotherByIndex with m_chargeConjugateDaughters=false, m_inclusiveDecay =fa...
Definition MCParticle.cc:615
ROOT::VecOps::RVec< float > get_vertex_y(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the production vertex y of the input MCParticles
Definition MCParticle.cc:242
ROOT::VecOps::RVec< float > get_pt(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the transverse momenta of the input MCParticles
Definition MCParticle.cc:176
std::vector< int > get_list_of_particles_from_decay(int i, ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
return the list of particles from the decay of a mother particle. i is the mother index in the Partic...
Definition MCParticle.cc:518
ROOT::VecOps::RVec< int > get_parentid(ROOT::VecOps::RVec< int > mcind, ROOT::VecOps::RVec< edm4hep::MCParticleData > mc, ROOT::VecOps::RVec< int > parents)
return the parent index of a given list of MC particles
Definition MCParticle.cc:431
ROOT::VecOps::RVec< float > get_e(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the energy of the input MCParticles
Definition MCParticle.cc:340
ROOT::VecOps::RVec< edm4hep::MCParticleData > mergeParticles(ROOT::VecOps::RVec< edm4hep::MCParticleData > x, ROOT::VecOps::RVec< edm4hep::MCParticleData > y)
concatenate both input vectors and return the resulting vector
Definition MCParticle.cc:184
ROOT::VecOps::RVec< float > get_time(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the time of the input MCParticles
Definition MCParticle.cc:194
ROOT::VecOps::RVec< float > AngleBetweenTwoMCParticles(ROOT::VecOps::RVec< edm4hep::MCParticleData > p1, ROOT::VecOps::RVec< edm4hep::MCParticleData > p2)
return the angle (3D) between two MCParticles :
Definition MCParticle.cc:674
ROOT::VecOps::RVec< float > get_vertex_x(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the production vertex x of the input MCParticles
Definition MCParticle.cc:234
ROOT::VecOps::RVec< float > get_endPoint_x(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the end point x of the input MCParticles
Definition MCParticle.cc:288
ROOT::VecOps::RVec< float > get_y(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the rapidity of the input MCParticles
Definition MCParticle.cc:392
ROOT::VecOps::RVec< float > get_p(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the momenta of the input MCParticles
Definition MCParticle.cc:350
std::vector< int > get_list_of_stable_particles_from_decay(int i, ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
return the list of stable particles from the decay of a mother particle, looking at the full decay ch...
Definition MCParticle.cc:484
std::vector< int > list_of_stable_particles_from_decay(int i, ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
obsolete: should use get_list_of_stable_particles_from_decay instead
Definition MCParticle.cc:548
ROOT::VecOps::RVec< float > get_pdg(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the PDG of the input MCParticles
Definition MCParticle.cc:202
int get_n(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the size of the input collection
Definition MCParticle.cc:422
ROOT::VecOps::RVec< edm4hep::Vector3d > get_endPoint(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the end point of the input MCParticles
Definition MCParticle.cc:258
ROOT::VecOps::RVec< float > get_mass(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the masses of the input MCParticles
Definition MCParticle.cc:312
ROOT::VecOps::RVec< float > get_charge(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the charges of the input MCParticles
Definition MCParticle.cc:384
ROOT::VecOps::RVec< float > get_vertex_z(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
return the production vertex z of the input MCParticles
Definition MCParticle.cc:250
FCC analyzers collection.
Definition Algorithms.h:15
Filter events based on a MCParticles PDGID.
Definition MCParticle.h:26
bool m_abs
Definition MCParticle.h:29
filter_pdgID(int arg_pdgid, bool arg_abs)
Definition MCParticle.cc:76
int m_pdgid
Definition MCParticle.h:28
bool operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
Definition MCParticle.cc:77
return the event primary vertex position and time (mm)
Definition MCParticle.h:79
TLorentzVector operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
Definition MCParticle.cc:103
get_EventPrimaryVertexP4()
Definition MCParticle.cc:102
return the event primary vertex (mm)
Definition MCParticle.h:72
get_EventPrimaryVertex(int arg_genstatus)
Definition MCParticle.cc:86
TVector3 operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
Definition MCParticle.cc:87
int m_genstatus
Definition MCParticle.h:74
get the decay of a given particle
Definition MCParticle.h:63
get_decay(int arg_mother, int arg_daughters, bool arg_inf)
Definition MCParticle.cc:42
bool operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
Definition MCParticle.cc:43
int m_mother
Definition MCParticle.h:65
int m_daughters
Definition MCParticle.h:66
bool m_inf
Definition MCParticle.h:67
A shorthand for get_indices, with m_chargeConjugateDaughters=false, inclusiveDecay=false.
Definition MCParticle.h:99
get_indices_ExclusiveDecay(int pdg_mother, std::vector< int > pdg_daughters, bool stableDaughters, bool chargeConjugate)
Definition MCParticle.cc:668
return a list of indices that correspond to a given MC decay. The list contains the index of the moth...
Definition MCParticle.h:87
get_indices(int pdg_mother, std::vector< int > pdg_daughters, bool stableDaughters, bool chargeConjugateMother, bool chargeConjugateDaughters, bool inclusiveDecay)
Definition MCParticle.cc:631
bool m_chargeConjugateMother
Definition MCParticle.h:92
bool m_chargeConjugateDaughters
Definition MCParticle.h:93
bool m_stableDaughters
Definition MCParticle.h:91
std::vector< int > m_pdg_daughters
Definition MCParticle.h:90
bool m_inclusiveDecay
Definition MCParticle.h:94
int m_pdg_mother
Definition MCParticle.h:89
ROOT::VecOps::RVec< int > operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
Definition MCParticle.cc:640
get MC history tree for a given MCParticle index
Definition MCParticle.h:56
get_tree(int arg_index)
Definition MCParticle.cc:133
int m_index
Definition MCParticle.h:58
ROOT::VecOps::RVec< int > operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in, ROOT::VecOps::RVec< int > ind)
Definition MCParticle.cc:134
select MCParticles with their status
Definition MCParticle.h:41
ROOT::VecOps::RVec< edm4hep::MCParticleData > operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
Definition MCParticle.cc:12
int m_status
Definition MCParticle.h:43
sel_genStatus(int arg_status)
Definition MCParticle.cc:11
select MCParticles with their PDG id
Definition MCParticle.h:48
int m_pdg
Definition MCParticle.h:50
ROOT::VecOps::RVec< edm4hep::MCParticleData > operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
Definition MCParticle.cc:25
bool m_chargeconjugate
Definition MCParticle.h:51
sel_pdgID(int arg_pdg, bool arg_chargeconjugate)
Definition MCParticle.cc:24
select MCParticles with transverse momentum greater than a minimum value [GeV]
Definition MCParticle.h:34
ROOT::VecOps::RVec< edm4hep::MCParticleData > operator()(ROOT::VecOps::RVec< edm4hep::MCParticleData > in)
Definition MCParticle.cc:63
sel_pt(float arg_min_pt)
Definition MCParticle.cc:62
float m_min_pt
Definition MCParticle.h:36