FCCAnalyses
Loading...
Searching...
No Matches
SmearObjects.h
Go to the documentation of this file.
1#ifndef SMEARING_ANALYZERS_H
2#define SMEARING_ANALYZERS_H
3
4#include <cmath>
5#include <vector>
6
8#include "ROOT/RVec.hxx"
9#include "TLorentzVector.h"
10#include "TMatrixDSym.h"
11#include "TRandom.h"
12#include "edm4hep/MCParticleData.h"
13#include <TMath.h>
14
15namespace FCCAnalyses {
16
17namespace SmearObjects {
18
21TVectorD TrackParamFromMC_DelphesConv(edm4hep::MCParticleData aMCParticle);
22
25 bool m_debug;
26 TRandom m_random;
28 SmearedTracks(float smear_d0, float smear_phi, float smear_omega,
29 float smear_z0, float smear_tlambda, bool debug);
30 ROOT::VecOps::RVec<edm4hep::TrackState>
31 operator()(const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
32 &allRecoParticles,
33 const ROOT::VecOps::RVec<edm4hep::TrackState> &alltracks,
34 const ROOT::VecOps::RVec<int> &RP2MC_indices,
35 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
36};
37
40ROOT::VecOps::RVec<edm4hep::TrackState> mcTrackParameters(
41 const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
42 &allRecoParticles,
43 const ROOT::VecOps::RVec<edm4hep::TrackState> &alltracks,
44 const ROOT::VecOps::RVec<int> &RP2MC_indices,
45 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
46
49TVectorD CovSmear(TVectorD x, TMatrixDSym C, TRandom *ran, bool debug);
50
54 bool m_debug;
55 TRandom m_random;
56 float m_scale;
57 SmearedTracksdNdx(float m_scale, bool debug);
58 ROOT::VecOps::RVec<edm4hep::Quantity>
59 operator()(const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
60 &allRecoParticles,
61 const ROOT::VecOps::RVec<edm4hep::Quantity> &dNdx,
62 const ROOT::VecOps::RVec<float> &length,
63 const ROOT::VecOps::RVec<int> &RP2MC_indices,
64 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
65};
66
69 bool m_debug;
70 TRandom m_random;
71 float m_scale;
72 SmearedTracksTOF(float m_scale, bool debug);
73 ROOT::VecOps::RVec<edm4hep::TrackerHit3DData>
74 operator()(const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
75 &allRecoParticles,
76 const ROOT::VecOps::RVec<edm4hep::TrackData> &trackdata,
77 const ROOT::VecOps::RVec<edm4hep::TrackerHit3DData> &trackerhits,
78 const ROOT::VecOps::RVec<float> &length,
79 const ROOT::VecOps::RVec<int> &RP2MC_indices,
80 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
81};
82
85 bool m_debug;
86 float m_scale;
87 int m_type;
88 int m_mode;
89 SmearedReconstructedParticle(float scale, int type, int mode, bool debug);
90
91 ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
92 operator()(const ROOT::VecOps::RVec<edm4hep::ReconstructedParticleData>
93 &allRecoParticles,
94 const ROOT::VecOps::RVec<int> &RP2MC_indices,
95 const ROOT::VecOps::RVec<edm4hep::MCParticleData> &mcParticles);
96};
97
98} // namespace SmearObjects
99} // namespace FCCAnalyses
100
101#endif
TVectorD TrackParamFromMC_DelphesConv(edm4hep::MCParticleData aMCParticle)
for a given MC particle, returns a "track state", i.e.
Definition SmearObjects.cc:22
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:185
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:240
FCC analyzers collection.
Definition Algorithms.h:15
generates new reco particles, smeared by given parameters
Definition SmearObjects.h:84
ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > operator()(const ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > &allRecoParticles, const ROOT::VecOps::RVec< int > &RP2MC_indices, const ROOT::VecOps::RVec< edm4hep::MCParticleData > &mcParticles)
Definition SmearObjects.cc:517
SmearedReconstructedParticle(float scale, int type, int mode, bool debug)
Definition SmearObjects.cc:497
generates new tracker hits, by rescaling the timing measurement
Definition SmearObjects.h:68
ROOT::VecOps::RVec< edm4hep::TrackerHit3DData > operator()(const ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > &allRecoParticles, const ROOT::VecOps::RVec< edm4hep::TrackData > &trackdata, const ROOT::VecOps::RVec< edm4hep::TrackerHit3DData > &trackerhits, const ROOT::VecOps::RVec< float > &length, const ROOT::VecOps::RVec< int > &RP2MC_indices, const ROOT::VecOps::RVec< edm4hep::MCParticleData > &mcParticles)
Definition SmearObjects.cc:400
bool m_debug
Definition SmearObjects.h:69
SmearedTracksTOF(float m_scale, bool debug)
Definition SmearObjects.cc:391
TRandom m_random
Definition SmearObjects.h:70
float m_scale
Definition SmearObjects.h:71
generates new track states, by rescaling the covariance matrix of the tracks
Definition SmearObjects.h:24
SmearedTracks(float smear_d0, float smear_phi, float smear_omega, float smear_z0, float smear_tlambda, bool debug)
Definition SmearObjects.cc:36
ROOT::VecOps::RVec< edm4hep::TrackState > operator()(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)
Definition SmearObjects.cc:49
bool m_debug
Definition SmearObjects.h:25
float m_smear_parameters[5]
Definition SmearObjects.h:27
TRandom m_random
Definition SmearObjects.h:26
generates new track dNdx, by rescaling the poisson error of the cluster count
Definition SmearObjects.h:53
ROOT::VecOps::RVec< edm4hep::Quantity > operator()(const ROOT::VecOps::RVec< edm4hep::ReconstructedParticleData > &allRecoParticles, const ROOT::VecOps::RVec< edm4hep::Quantity > &dNdx, const ROOT::VecOps::RVec< float > &length, const ROOT::VecOps::RVec< int > &RP2MC_indices, const ROOT::VecOps::RVec< edm4hep::MCParticleData > &mcParticles)
Definition SmearObjects.cc:313
bool m_debug
Definition SmearObjects.h:54
SmearedTracksdNdx(float m_scale, bool debug)
Definition SmearObjects.cc:304
TRandom m_random
Definition SmearObjects.h:55
float m_scale
Definition SmearObjects.h:56