1#ifndef ANALYZERS_SOURCE_RECONSTRUCTED_PARTICLE_H
2#define ANALYZERS_SOURCE_RECONSTRUCTED_PARTICLE_H
5#include "ROOT/RVec.hxx"
8#include "edm4hep/RecoMCParticleLinkCollection.h"
9#include "edm4hep/ReconstructedParticleCollection.h"
25 explicit selPDG(
const int pdgID);
33 edm4hep::ReconstructedParticleCollection
34 operator()(
const edm4hep::RecoMCParticleLinkCollection &inLinkColl);
57 edm4hep::ReconstructedParticleCollection
58 operator()(
const edm4hep::RecoMCParticleLinkCollection &inLinkColl);
72 explicit selPt(
float minPt);
79 edm4hep::ReconstructedParticleCollection
80 operator()(
const edm4hep::ReconstructedParticleCollection &inColl);
93 explicit selUpTo(
const size_t size);
100 edm4hep::ReconstructedParticleCollection
101 operator()(
const edm4hep::ReconstructedParticleCollection &inColl);
122 edm4hep::ReconstructedParticleCollection
123 operator()(
const edm4hep::RecoMCParticleLinkCollection &inLinkColl);
134ROOT::VecOps::RVec<float>
135getP(
const edm4hep::ReconstructedParticleCollection &inColl);
143ROOT::VecOps::RVec<float>
144getPt(
const edm4hep::ReconstructedParticleCollection &inColl);
152ROOT::VecOps::RVec<float>
153getY(
const edm4hep::ReconstructedParticleCollection &inColl);
161ROOT::VecOps::RVec<float>
162getE(
const edm4hep::ReconstructedParticleCollection &inColl);
170ROOT::VecOps::RVec<float>
171getMass(
const edm4hep::ReconstructedParticleCollection &inColl);
179ROOT::VecOps::RVec<float>
180getCharge(
const edm4hep::ReconstructedParticleCollection &inColl);
190edm4hep::ReconstructedParticleCollection
191sortByPt(
const edm4hep::ReconstructedParticleCollection &inColl);
206edm4hep::ReconstructedParticleCollection
207remove(
const edm4hep::ReconstructedParticleCollection &inColl,
208 const edm4hep::ReconstructedParticle &inPartToBeRemoved,
209 const bool matching =
false);
225edm4hep::ReconstructedParticleCollection
226remove(
const edm4hep::ReconstructedParticleCollection &inColl,
227 const edm4hep::ReconstructedParticleCollection &inPartsToBeRemoved,
228 const bool matching =
false);
237edm4hep::ReconstructedParticleCollection
238merge(
const edm4hep::ReconstructedParticleCollection &inColl1,
239 const edm4hep::ReconstructedParticleCollection &inColl2);
260 edm4hep::ReconstructedParticleCollection
261 operator()(
const edm4hep::ReconstructedParticleCollection &inColl);
282 edm4hep::ReconstructedParticleCollection
283 operator()(
const edm4hep::ReconstructedParticleCollection &inColl);
edm4hep::ReconstructedParticleCollection sortByPt(const edm4hep::ReconstructedParticleCollection &inColl)
Sort input particles by pT.
Definition ReconstructedParticleSource.cc:195
ROOT::VecOps::RVec< float > getPt(const edm4hep::ReconstructedParticleCollection &inColl)
Get transverse momenta (pT) of the input particles.
Definition ReconstructedParticleSource.cc:134
ROOT::VecOps::RVec< float > getP(const edm4hep::ReconstructedParticleCollection &inColl)
Get momenta of the input reconstructed particles.
Definition ReconstructedParticleSource.cc:122
ROOT::VecOps::RVec< float > getY(const edm4hep::ReconstructedParticleCollection &inColl)
Get rapidity (y) of the input reconstructed particles.
Definition ReconstructedParticleSource.cc:147
ROOT::VecOps::RVec< float > getE(const edm4hep::ReconstructedParticleCollection &inColl)
Get energy (E) of the input reconstructed particles.
Definition ReconstructedParticleSource.cc:160
ROOT::VecOps::RVec< float > getCharge(const edm4hep::ReconstructedParticleCollection &inColl)
Get charge of the input reconstructed particles.
Definition ReconstructedParticleSource.cc:183
FCC analyzers collection.
Definition Algorithms.h:15
Build the recoil from an arbitrary list of input resonances at the specified center of mass energy.
Definition ReconstructedParticleSource.h:268
const float m_sqrts
Definition ReconstructedParticleSource.h:269
Build two particle resonances from an arbitrary list of input reconstructed particles.
Definition ReconstructedParticleSource.h:246
const float m_resonanceMass
Definition ReconstructedParticleSource.h:247
Analyzer to select reconstructed particles associated with MC particle with the specified absolute va...
Definition ReconstructedParticleSource.h:41
const int m_absPdg
Definition ReconstructedParticleSource.h:42
Analyzer to select reconstructed particles associated with the MC particle of the desired generator s...
Definition ReconstructedParticleSource.h:108
const int m_status
Definition ReconstructedParticleSource.h:109
Analyzer to select reconstructed particles associated with MC particle of the specified PDG ID.
Definition ReconstructedParticleSource.h:18
const int m_pdg
Definition ReconstructedParticleSource.h:19
Select reconstructed particles with transverse momentum greater than a minimum value [GeV].
Definition ReconstructedParticleSource.h:65
const float m_minPt
Definition ReconstructedParticleSource.h:66
Analyzer to select specified number of reconstructed particles.
Definition ReconstructedParticleSource.h:86
const size_t m_size
Definition ReconstructedParticleSource.h:87