1#ifndef SMEAR_OBJECTS_ANALYZERS_H
2#define SMEAR_OBJECTS_ANALYZERS_H
8#include "ROOT/RVec.hxx"
9#include "TLorentzVector.h"
11#include "TMatrixDSym.h"
14#include "edm4hep/MCParticleData.h"
15#include "edm4hep/RecDqdxData.h"
30 float m_smear_parameters[5];
31 SmearedTracks(
float smear_d0,
float smear_phi,
float smear_omega,
32 float smear_z0,
float smear_tlambda,
bool debug);
33 ROOT::VecOps::RVec<edm4hep::TrackState>
34 operator()(
const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
36 const ROOT::VecOps::RVec<edm4hep::TrackState> &alltracks,
37 const ROOT::VecOps::RVec<int> &RP2MC_indices,
38 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
44 const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
46 const ROOT::VecOps::RVec<edm4hep::TrackState> &alltracks,
47 const ROOT::VecOps::RVec<int> &RP2MC_indices,
48 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
52TVectorD
CovSmear(TVectorD x, TMatrixDSym C, TRandom *ran,
bool debug);
74 ROOT::VecOps::RVec<edm4hep::RecDqdxData>
75 operator()(
const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
77 const TrackUtils::TrackDqdxHandler &dNdxHandler,
78 const ROOT::VecOps::RVec<float> &length,
79 const ROOT::VecOps::RVec<int> &RP2MC_indices,
80 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
89 ROOT::VecOps::RVec<edm4hep::TrackerHit3DData>
90 operator()(
const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
92 const ROOT::VecOps::RVec<edm4hep::TrackData> &trackdata,
93 const ROOT::VecOps::RVec<edm4hep::TrackerHit3DData> &trackerhits,
94 const ROOT::VecOps::RVec<float> &length,
95 const ROOT::VecOps::RVec<int> &RP2MC_indices,
96 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
107 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
108 operator()(
const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
110 const ROOT::VecOps::RVec<int> &RP2MC_indices,
111 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
TVectorD TrackParamFromMC_DelphesConv(edm4hep::MCParticleData aMCParticle)
for a given MC particle, returns a "track state", i.e.
Definition SmearObjects.cc:20
ROOT::VecOps::RVec< edm4hep::TrackState > mcTrackParameters(const ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > &allRecoParticles, const ROOT::VecOps::RVec< edm4hep::TrackState > &alltracks, const ROOT::VecOps::RVec< int > &RP2MC_indices, const ROOT::VecOps::RVec< edm4hep::MCParticleData > &mcParticles)
used to validate the method above.
Definition SmearObjects.cc:183
TVectorD CovSmear(TVectorD x, TMatrixDSym C, TRandom *ran, bool debug)
generates random values for a vector, given the covariance matrix of its components,...
Definition SmearObjects.cc:238
FCC analyzers collection.
Definition Algorithms.h:15
generates new reco particles, smeared by given parameters
Definition SmearObjects.h:100
bool m_debug
Definition SmearObjects.h:101
int m_type
Definition SmearObjects.h:103
float m_scale
Definition SmearObjects.h:102
int m_mode
Definition SmearObjects.h:104
generates new tracker hits, by rescaling the timing measurement
Definition SmearObjects.h:84
TRandom m_random
Definition SmearObjects.h:86
bool m_debug
Definition SmearObjects.h:85
float m_scale
Definition SmearObjects.h:87
generates new track states, by rescaling the covariance matrix of the tracks
Definition SmearObjects.h:27
TRandom m_random
Definition SmearObjects.h:29
bool m_debug
Definition SmearObjects.h:28
Generates new track dNdx, by rescaling the Poisson error of the cluster count.
Definition SmearObjects.h:58
TRandom m_random
Debug flag.
Definition SmearObjects.h:60
bool m_debug
Definition SmearObjects.h:59
float m_scale
Definition SmearObjects.h:61