EIC Software
Reference for
EIC
simulation and reconstruction software on GitHub
Home page
Related Pages
Modules
Namespaces
Classes
Files
External Links
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
SmearTrackingPreview_0_2_B3T.cxx
Go to the documentation of this file.
Or view
the newest version in sPHENIX GitHub for file SmearTrackingPreview_0_2_B3T.cxx
1
// Based on
2
// https://physdiv.jlab.org/DetectorMatrix/
3
// Using everything from v 0.1 from June 16 2020 except:
4
// - EMCal coverage reduced to |eta|<3.5
5
// - Tracking parameters from B=3T Matrix preview at
6
// https://indico.bnl.gov/event/9984/contributions/43066/attachments/31173/49186/YR_Detector_Matrix_Tracking_only_10282020.xlsx
7
8
// Reminder:
9
// Acceptance::Zone(double theta = 0., double = TMath::Pi(),
10
// double phi = 0., double = TMath::TwoPi(),
11
// double E = 0., double = TMath::Infinity(),
12
// double p = 0., double = TMath::Infinity(),
13
// double pt = 0., double = TMath::Infinity(),
14
// double pz = -TMath::Infinity(), double = TMath::Infinity());
15
16
17
#include "
eicsmear/erhic/VirtualParticle.h
"
18
#include "
eicsmear/smear/Acceptance.h
"
19
#include "
eicsmear/smear/Device.h
"
20
#include "
eicsmear/smear/Detector.h
"
21
#include "
eicsmear/smear/Smearer.h
"
22
#include "
eicsmear/smear/ParticleMCS.h
"
23
#include "
eicsmear/smear/PerfectID.h
"
24
#include <
eicsmear/smear/Smear.h
>
25
#include <
eicsmear/erhic/ParticleMC.h
>
26
#include "Math/Vector4D.h"
27
28
// declare static --> local to this file, won't clash with others
29
static
double
ThetaFromEta
(
const
double
eta
);
30
31
Smear::Detector
BuildTrackingPreview_0_2_B3T
() {
32
gSystem->Load(
"libeicsmear"
);
33
34
// Create the detector object to hold all devices
35
Smear::Detector
det;
36
37
// The framework provides implementations of three kinematic calculation methods
38
// from smeared values
39
// NM - "Null method" uses the scattered lepton.
40
// DA - Double Angle method
41
// JB - Jacquet-Blondel method
42
// Important Notes:
43
// - All methods rely on measured energy and momentum components (lepton for NM, hadrons for DA, JB).
44
// In order to rely on these methods, you therefore _need_ to cover the relevant
45
// regions with smearers for momentum, angles and energy.
46
// - This is exacerbated by the fact that in the current implementation a smearing of e.g. P
47
// may fool the framework into assuming theta is measured to be the initialization value 0,
48
// leading to nan or inf or otherwise wrong calculations.
49
// NM catches P=0 and replaces it in calculations with E,
50
// but JB and DA either don't work at all or report nonsenical values.
51
// - It may be advantageous to use a P measurement in place of E in a more general case.
52
// In the future, we plan to change to a weighted mean approach.
53
// The summary of the above is that the user should be mindful of the limitations and assumptions in
54
// the calculation of smeared kinematics. Sophisticated calculations are best redone in user code
55
// from the smeared particles directly.
56
det.
SetEventKinematicsCalculator
(
"NM DA JB"
);
57
58
// Perfect phi and theta for all particles
59
// ---------------------------------------
60
// TODO: Better options?
61
// total coverage of the handbook for tracker, ecal, and hcal is -3.5 < eta < 3.5
62
// muons no longer need special treatment! They are part of kAll
63
Smear::Acceptance::Zone
AngleZoneCommon(
ThetaFromEta
( 3.5 ),
ThetaFromEta
( -3.5 ));
64
Smear::Device
SmearThetaCommon(
Smear::kTheta
,
"0.0"
);
65
SmearThetaCommon.
Accept
.
AddZone
(AngleZoneCommon);
66
SmearThetaCommon.
Accept
.
SetGenre
(
Smear::kAll
);
67
det.
AddDevice
(SmearThetaCommon);
68
69
Smear::Device
SmearPhiCommon(
Smear::kPhi
,
"0.0"
);
70
SmearPhiCommon.
Accept
.
AddZone
(AngleZoneCommon);
71
SmearPhiCommon.
Accept
.
SetGenre
(
Smear::kAll
);
72
det.
AddDevice
(SmearPhiCommon);
73
74
75
// Tracking (B = 3 T)
76
// --------
77
// Note: Smear::kCharged checks pdg charge, so includes muons (good)
78
// eta = -3.5 -- -2.5
79
// sigma_p/p ~ 0.1% p + 2%
80
Smear::Acceptance::Zone
TrackBack1Zone(
ThetaFromEta
( -2.5 ),
ThetaFromEta
( -3.5 ));
81
Smear::Device
TrackBack1P(
Smear::kP
,
"sqrt( pow ( 0.001*P*P, 2) + pow ( 0.02*P, 2) )"
);
82
TrackBack1P.
Accept
.
AddZone
(TrackBack1Zone);
83
TrackBack1P.
Accept
.
SetCharge
(
Smear::kCharged
);
84
det.
AddDevice
(TrackBack1P);
85
86
// eta = -2.5 -- -1
87
// sigma_p/p ~ 0.02% p + 1%
88
Smear::Acceptance::Zone
TrackBack2Zone(
ThetaFromEta
( -1 ),
ThetaFromEta
( -2.5 ));
89
Smear::Device
TrackBack2P(
Smear::kP
,
"sqrt( pow ( 0.0002*P*P, 2) + pow ( 0.01*P, 2) )"
);
90
TrackBack2P.
Accept
.
AddZone
(TrackBack2Zone);
91
TrackBack2P.
Accept
.
SetCharge
(
Smear::kCharged
);
92
det.
AddDevice
(TrackBack2P);
93
94
// TODO: What to make of this?
95
// "200 MeV/c with 50% acceptance (similar for pi and K)"
96
// Add efficiency?
97
// eta = -1 -- +1
98
// sigma_p/p ~ 0.02% p + 0.05%
99
Smear::Acceptance::Zone
TrackBarrelZone(
ThetaFromEta
( 1 ),
ThetaFromEta
( -1 ));
100
Smear::Device
TrackBarrelP(
Smear::kP
,
"sqrt( pow ( 0.0002*P*P, 2) + pow ( 0.005*P, 2) )"
);
101
TrackBarrelP.
Accept
.
AddZone
(TrackBarrelZone);
102
TrackBarrelP.
Accept
.
SetCharge
(
Smear::kCharged
);
103
det.
AddDevice
(TrackBarrelP);
104
105
// eta = 1 -- 2.5
106
// sigma_p/p ~ 0.02% p+1%
107
Smear::Acceptance::Zone
TrackFwd2Zone(
ThetaFromEta
( 2.5 ),
ThetaFromEta
( 1 ));
108
Smear::Device
TrackFwd2P(
Smear::kP
,
"sqrt( pow ( 0.0002*P*P, 2) + pow ( 0.01*P, 2) )"
);
109
TrackFwd2P.
Accept
.
AddZone
(TrackFwd2Zone);
110
TrackFwd2P.
Accept
.
SetCharge
(
Smear::kCharged
);
111
det.
AddDevice
(TrackFwd2P);
112
113
// eta = 2.5 -- 3.5
114
// sigma_p/p ~ 0.1% p+2%
115
Smear::Acceptance::Zone
TrackFwd1Zone(
ThetaFromEta
( 3.5 ),
ThetaFromEta
( 2.5 ));
116
Smear::Device
TrackFwd1P(
Smear::kP
,
"sqrt( pow ( 0.001*P*P, 2) + pow ( 0.02*P, 2) )"
);
117
TrackFwd1P.
Accept
.
AddZone
(TrackFwd1Zone);
118
TrackFwd1P.
Accept
.
SetCharge
(
Smear::kCharged
);
119
det.
AddDevice
(TrackFwd1P);
120
121
// EM Calorimeters
122
// ---------------
123
// Note: Smear::kElectromagnetic == gamma + e. Does not include muons (good)
124
125
// Calorimeter resolution usually given as sigma_E/E = const% + stocastic%/Sqrt{E}
126
// EIC Smear needs absolute sigma: sigma_E = Sqrt{const*const*E*E + stoc*stoc*E}
127
128
// Back
129
// eta = -3.5 -- -2
130
// stoch. = 2%
131
Smear::Acceptance::Zone
EmcalBackZone(
ThetaFromEta
( -2 ),
ThetaFromEta
( -3.5 ));
132
Smear::Device
EmcalBack(
Smear::kE
,
"sqrt( pow ( 0.0*E,2 ) + pow( 0.02,2)*E)"
);
133
EmcalBack.
Accept
.
AddZone
(EmcalBackZone);
134
EmcalBack.
Accept
.
SetGenre
(
Smear::kElectromagnetic
);
135
det.
AddDevice
(EmcalBack);
136
137
// MidBack
138
// eta = -2 -- -1
139
// stoch. = 7%
140
Smear::Acceptance::Zone
EmcalMidBackZone(
ThetaFromEta
( -1 ),
ThetaFromEta
( -2 ));
141
Smear::Device
EmcalMidBack(
Smear::kE
,
"sqrt( pow ( 0.0*E,2 ) + pow( 0.07,2)*E)"
);
142
EmcalMidBack.
Accept
.
AddZone
(EmcalMidBackZone);
143
EmcalMidBack.
Accept
.
SetGenre
(
Smear::kElectromagnetic
);
144
det.
AddDevice
(EmcalMidBack);
145
146
// Forward
147
// eta = -1 -- 3.5
148
// stoch. = 10-12%, use 12%
149
Smear::Acceptance::Zone
EmcalFwdZone(
ThetaFromEta
( 3.5 ),
ThetaFromEta
( -1 ));
150
Smear::Device
EmcalFwd(
Smear::kE
,
"sqrt( pow ( 0.0*E,2 ) + pow( 0.12,2)*E)"
);
151
EmcalFwd.
Accept
.
AddZone
(EmcalFwdZone);
152
EmcalFwd.
Accept
.
SetGenre
(
Smear::kElectromagnetic
);
153
det.
AddDevice
(EmcalFwd);
154
155
// Could turn on perfect PID
156
// Make sure to not cover more than is covered by the other detectors.
157
// Using the smallest region here, but we could add a second one for
158
// the extended EmCal range
159
// Smear::Acceptance::Zone acceptpid(ThetaFromEta(3.5),ThetaFromEta(-3.5));
160
// Smear::PerfectID pid;
161
// pid.Accept.AddZone(acceptpid);
162
// det.AddDevice( pid );
163
164
// Hadronic Calorimeters
165
// ----------------------
166
// Note: kHadronic == |pdg|>110.
167
168
// Back
169
// eta = -3.5 -- -1
170
// stoch. = 50%
171
Smear::Acceptance::Zone
HcalBackZone(
ThetaFromEta
( -1 ),
ThetaFromEta
( -3.5 ));
172
Smear::Device
HcalBack(
Smear::kE
,
"sqrt(pow( 0.0*E, 2) + pow ( 0.5,2) *E)"
);
173
HcalBack.
Accept
.
AddZone
(HcalBackZone);
174
HcalBack.
Accept
.
SetGenre
(
Smear::kHadronic
);
175
det.
AddDevice
(HcalBack);
176
177
// Barrel
178
// eta = -1 -- 1
179
// The matrix has nothing. As examples, one could turn to
180
// ~CMS
181
// Smear::Device HcalBarrel(Smear::kE, "sqrt( pow( 0.07*E, 2) + pow( 0.85, 2)*E)");
182
// ~Zeus
183
// Smear::Device HcalBarrel(Smear::kE, "sqrt( pow( 0.02*E, 2) + pow( 0.35,2) *E)");
184
185
// Smear::Acceptance::Zone HcalBarrelZone(ThetaFromEta ( 1 ),ThetaFromEta ( -1 ));
186
// HcalBarrel.Accept.AddZone(HcalBarrelZone);
187
// HcalBarrel.Accept.SetGenre(Smear::kHadronic);
188
// det.AddDevice(HcalBarrel);
189
190
// Forward
191
// eta = 1 -- 3.5
192
// stoch. = 50%
193
Smear::Acceptance::Zone
HcalFwdZone(
ThetaFromEta
( 3.5 ),
ThetaFromEta
( 1 ));
194
Smear::Device
HcalFwd(
Smear::kE
,
"sqrt(pow( 0.0*E, 2) + pow ( 0.5,2) *E)"
);
195
HcalFwd.
Accept
.
AddZone
(HcalFwdZone);
196
HcalFwd.
Accept
.
SetGenre
(
Smear::kHadronic
);
197
det.
AddDevice
(HcalFwd);
198
199
// Not covered:
200
// Tracker material budget X/X0 <~5%
201
// Low-Q^2 tagger: -6.9<eta<-5.8: Delta_theta/theta < 1.5%; 10^-6 < Q2 < 10^-2 GeV2
202
// Proton spectrometer: eta>6.2: sigma_intrinsic(|t|)/|t| < 1%; Acceptance: 0.2 < pT < 1.2 GeV/c
203
// Barrel vertexing: sigma_xyz ~ 20 microns, d0(z) ~ d0(r phi) ~ (20 microns)/(pT [GeV]) + 5 microns
204
205
return
det;
206
}
207
208
// -------------------------------------------------------------------
209
double
ThetaFromEta
(
const
double
eta
) {
210
if
( !std::isnan(eta) && !std::isinf(eta) ) {
211
return
2.0 * atan( exp( -eta ));
212
}
213
throw
std::runtime_error(
"ThetaFromEta called with NaN or Inf"
);
214
return
-1;
215
}
eicsmeardetectors
blob
master
SmearTrackingPreview_0_2_B3T.cxx
Built by
Jin Huang
. updated:
Mon Jan 22 2024 12:43:37
using
1.8.2 with
EIC GitHub integration