EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SmearMatrixDetector_0_1_TOF.cxx
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file SmearMatrixDetector_0_1_TOF.cxx
1 // Based on
2 // https://physdiv.jlab.org/DetectorMatrix/
3 // From June 16 2020
4 
5 // W.r.t. to the [Detector Requirements and R&D Handbook](http://www.eicug.org/web/sites/default/files/EIC_HANDBOOK_v1.2.pdf) there have been three changes.
6 
7 // For the Backward Detector the Tracking Resolution column was updated as follows:
8 // eta=-3.5 - -2.5: sigma_p/p ~ 0.1%×p+2.0% -> sigma_p/p ~ 0.1%×p ⊕ 0.5%
9 // eta=-2.5 - -2.0: sigma_p/p ~ 0.1%×p+1.0% -> sigma_p/p ~ 0.1%×p ⊕ 0.5%
10 // eta=-2.0 - -1.0: sigma_p/p ~ 0.1%×p+1.0% -> sigma_p/p ~ 0.05%×p ⊕ 0.5%
11 // The abstract, reference files, and note sections of these fields have been updated accordingly.
12 
13 // Important Notes:
14 // - The implementation of the original handbook matrix in SmearHandBook_1_2.cxx
15 // took some liberties and added constant terms and other specifications where
16 // they were/are missing in the matrix.
17 // The implementation here explicitly does not do so, and aims to be completely faithful
18 // to the configuration available online.
19 // - Where ranges are given, the more conservative number is chosen.
20 // - Without available specifications, angular resolution is assumed to be perfect.
21 // -
22 
25 #include "eicsmear/smear/Device.h"
27 #include "eicsmear/smear/Smearer.h"
30 #include <eicsmear/smear/Smear.h>
33 #include "Math/Vector4D.h"
34 
35 // declare static --> local to this file, won't clash with others
36 static double ThetaFromEta( const double eta );
37 
39  gSystem->Load("libeicsmear");
40 
41  // Create the detector object to hold all devices
42  Smear::Detector det;
43 
44  // The framework provides implementations of three kinematic calculation methods
45  // from smeared values
46  // NM - "Null method" uses the scattered lepton.
47  // DA - Double Angle method
48  // JB - Jacquet-Blondel method
49  // Important Notes:
50  // - All methods rely on measured energy and momentum components (lepton for NM, hadrons for DA, JB).
51  // In order to rely on these methods, you therefore _need_ to cover the relevant
52  // regions with smearers for momentum, angles and energy.
53  // - This is exacerbated by the fact that in the current implementation a smearing of e.g. P
54  // may fool the framework into assuming theta is measured to be the initialization value 0,
55  // leading to nan or inf or otherwise wrong calculations.
56  // NM catches P=0 and replaces it in calculations with E,
57  // but JB and DA either don't work at all or report nonsenical values.
58  // - It may be advantageous to use a P measurement in place of E in a more general case.
59  // In the future, we plan to change to a weighted mean approach.
60  // The summary of the above is that the user should be mindful of the limitations and assumptions in
61  // the calculation of smeared kinematics. Sophisticated calculations are best redone in user code
62  // from the smeared particles directly.
63  det.SetEventKinematicsCalculator("NM DA JB");
64 
65  // IMPORTANT: There are two traps (will be addressed in future releases):
66  // 1) If you smear E but don't provide phi, theta smearers, those values will be
67  // set to 0, not to a fault value and not to the truth level
68  // 2) If you do provide such a smearer, pt and pz will be changed
69  // by a consistency enforcer in Detector::Smear()
70 
71 
72  // Perfect phi and theta for all particles
73  // ---------------------------------------
74  // total coverage of the handbook for tracker and hcal is -3.5 < eta < 3.5
75  Smear::Acceptance::Zone AngleZoneHadronic(ThetaFromEta ( 3.5 ),ThetaFromEta ( -3.5 ));
76  Smear::Device SmearThetaHadronic(Smear::kTheta, "0.0");
77  SmearThetaHadronic.Accept.AddZone(AngleZoneHadronic);
78  SmearThetaHadronic.Accept.SetGenre(Smear::kHadronic);
79  det.AddDevice(SmearThetaHadronic);
80 
81  Smear::Device SmearPhiHadronic(Smear::kPhi, "0.0");
82  SmearPhiHadronic.Accept.AddZone(AngleZoneHadronic);
83  SmearPhiHadronic.Accept.SetGenre(Smear::kHadronic);
84  det.AddDevice(SmearPhiHadronic);
85 
86  // muons are neither hadrons nor electromgnetic
87  Smear::Acceptance::Zone AngleZoneMuon(ThetaFromEta ( 3.5 ),ThetaFromEta ( -3.5 ));
88  Smear::Device SmearThetaMuon(Smear::kTheta, "0.0");
89  SmearThetaMuon.Accept.AddZone(AngleZoneMuon);
90  SmearThetaMuon.Accept.AddParticle(13);
91  SmearThetaMuon.Accept.AddParticle(-13);
92  det.AddDevice(SmearThetaMuon);
93 
94  Smear::Device SmearPhiMuon(Smear::kPhi, "0.0");
95  SmearPhiMuon.Accept.AddZone(AngleZoneMuon);
96  SmearPhiMuon.Accept.AddParticle(13);
97  SmearPhiMuon.Accept.AddParticle(-13);
98  det.AddDevice(SmearPhiMuon);
99 
100  // emcal stretches to -4.5 < eta < 4.5
101  Smear::Acceptance::Zone AngleZoneEmcal(ThetaFromEta ( 4.5 ),ThetaFromEta ( -4.5 ));
102  Smear::Device SmearThetaEmcal(Smear::kTheta, "0.0");
103  SmearThetaEmcal.Accept.AddZone(AngleZoneEmcal);
104  SmearThetaEmcal.Accept.SetGenre(Smear::kElectromagnetic);
105  det.AddDevice(SmearThetaEmcal);
106 
107  Smear::Device SmearPhiEmcal(Smear::kPhi, "0.0");
108  SmearPhiEmcal.Accept.AddZone(AngleZoneEmcal);
109  SmearPhiEmcal.Accept.SetGenre(Smear::kElectromagnetic);
110  det.AddDevice(SmearPhiEmcal);
111 
112  // Tracking
113  // --------
114  // Note: Smear::kCharged checks pdg charge, so includes muons (good)
115  // eta = -3.5 -- -2.0
116  // sigma_p/p ~ 0.1% p+0.5%
117  Smear::Acceptance::Zone TrackBack1Zone(ThetaFromEta ( -2.5 ),ThetaFromEta ( -3.5 ));
118  Smear::Device TrackBack1P(Smear::kP, "sqrt( pow ( 0.001*P*P, 2) + pow ( 0.005*P, 2) )");
119  TrackBack1P.Accept.AddZone(TrackBack1Zone);
120  TrackBack1P.Accept.SetCharge(Smear::kCharged);
121  det.AddDevice(TrackBack1P);
122 
123  // eta = -2.0 -- -1
124  // sigma_p/p ~ 0.05% p+ 0.5%
125  Smear::Acceptance::Zone TrackBack2Zone(ThetaFromEta ( -1 ),ThetaFromEta ( -2.5 ));
126  Smear::Device TrackBack2P(Smear::kP, "sqrt( pow ( 0.0005*P*P, 2) + pow ( 0.005*P, 2) )");
127  TrackBack2P.Accept.AddZone(TrackBack2Zone);
128  TrackBack2P.Accept.SetCharge(Smear::kCharged);
129  det.AddDevice(TrackBack2P);
130 
131  // eta = -1 -- +1
132  // sigma_p/p ~ 0.05% p+0.5%
133  Smear::Acceptance::Zone TrackBarrelZone(ThetaFromEta ( 1 ),ThetaFromEta ( -1 ));
134  Smear::Device TrackBarrelP(Smear::kP, "sqrt( pow ( 0.0005*P*P, 2) + pow ( 0.005*P, 2) )");
135  TrackBarrelP.Accept.AddZone(TrackBarrelZone);
136  TrackBarrelP.Accept.SetCharge(Smear::kCharged);
137  det.AddDevice(TrackBarrelP);
138 
139  // eta = 1 -- 2.5
140  // sigma_p/p ~ 0.05% p+1.0%
141  Smear::Acceptance::Zone TrackFwd2Zone(ThetaFromEta ( 2.5 ),ThetaFromEta ( 1 ));
142  Smear::Device TrackFwd2P(Smear::kP, "sqrt( pow ( 0.0005*P*P, 2) + pow ( 0.01*P, 2) )");
143  TrackFwd2P.Accept.AddZone(TrackFwd2Zone);
144  TrackFwd2P.Accept.SetCharge(Smear::kCharged);
145  det.AddDevice(TrackFwd2P);
146 
147  // eta = 2.5 -- 3.5
148  // sigma_p/p ~ 0.1% p+2.0%
149  Smear::Acceptance::Zone TrackFwd1Zone(ThetaFromEta ( 3.5 ),ThetaFromEta ( 2.5 ));
150  Smear::Device TrackFwd1P(Smear::kP, "sqrt( pow ( 0.001*P*P, 2) + pow ( 0.02*P, 2) )");
151  TrackFwd1P.Accept.AddZone(TrackFwd1Zone);
152  TrackFwd1P.Accept.SetCharge(Smear::kCharged);
153  det.AddDevice(TrackFwd1P);
154 
155  // EM Calorimeters
156  // ---------------
157  // Note: Smear::kElectromagnetic == gamma + e. Does not include muons (good)
158 
159  // Calorimeter resolution usually given as sigma_E/E = const% + stocastic%/Sqrt{E}
160  // EIC Smear needs absolute sigma: sigma_E = Sqrt{const*const*E*E + stoc*stoc*E}
161 
162  // Back
163  // eta = -4.5 -- -2
164  // stoch. = 2%
165  Smear::Acceptance::Zone EmcalBackZone(ThetaFromEta ( -2 ),ThetaFromEta ( -4.5 ));
166  Smear::Device EmcalBack(Smear::kE, "sqrt( pow ( 0.0*E,2 ) + pow( 0.02,2)*E)");
167  EmcalBack.Accept.AddZone(EmcalBackZone);
169  det.AddDevice(EmcalBack);
170 
171  // MidBack
172  // eta = -2 -- -1
173  // stoch. = 7%
174  Smear::Acceptance::Zone EmcalMidBackZone(ThetaFromEta ( -1 ),ThetaFromEta ( -2 ));
175  Smear::Device EmcalMidBack(Smear::kE, "sqrt( pow ( 0.0*E,2 ) + pow( 0.07,2)*E)");
176  EmcalMidBack.Accept.AddZone(EmcalMidBackZone);
178  det.AddDevice(EmcalMidBack);
179 
180  // Forward
181  // eta = -1 -- 4.5
182  // stoch. = 10-12%, use 12%
183  Smear::Acceptance::Zone EmcalFwdZone(ThetaFromEta ( 4.5 ),ThetaFromEta ( -1 ));
184  Smear::Device EmcalFwd(Smear::kE, "sqrt( pow ( 0.0*E,2 ) + pow( 0.12,2)*E)");
185  EmcalFwd.Accept.AddZone(EmcalFwdZone);
187  det.AddDevice(EmcalFwd);
188 
189  // TODO: Add PID
190  // Could turn on perfect PID
191  // Make sure to not cover more than is covered by the other detectors.
192  // Using the smallest region here, but we could add a second one for
193  // the extended EmCal range
194  // Smear::Acceptance::Zone acceptpid(ThetaFromEta(3.5),ThetaFromEta(-3.5));
195  // Smear::PerfectID pid;
196  // pid.Accept.AddZone(acceptpid);
197  // det.AddDevice( pid );
198 
199  // Hadronic Calorimeters
200  // ----------------------
201  // Note: kHadronic == |pdg|>110.
202 
203  // Back
204  // eta = -3.5 -- -1
205  // stoch. = 50%
206  Smear::Acceptance::Zone HcalBackZone(ThetaFromEta ( -1 ),ThetaFromEta ( -3.5 ));
207  Smear::Device HcalBack(Smear::kE, "sqrt(pow( 0.0*E, 2) + pow ( 0.5,2) *E)");
208  HcalBack.Accept.AddZone(HcalBackZone);
209  HcalBack.Accept.SetGenre(Smear::kHadronic);
210  det.AddDevice(HcalBack);
211 
212  // Barrel
213  // eta = -1 -- 1
214  // The matrix has nothing. As examples, one could turn to
215  // ~CMS
216  // Smear::Device HcalBarrel(Smear::kE, "sqrt( pow( 0.07*E, 2) + pow( 0.85, 2)*E)");
217  // ~Zeus
218  // Smear::Device HcalBarrel(Smear::kE, "sqrt( pow( 0.02*E, 2) + pow( 0.35,2) *E)");
219 
220  // Smear::Acceptance::Zone HcalBarrelZone(ThetaFromEta ( 1 ),ThetaFromEta ( -1 ));
221  // HcalBarrel.Accept.AddZone(HcalBarrelZone);
222  // HcalBarrel.Accept.SetGenre(Smear::kHadronic);
223  // det.AddDevice(HcalBarrel);
224 
225  // Forward
226  // eta = 1 -- 3.5
227  // stoch. = 50%
228  Smear::Acceptance::Zone HcalFwdZone(ThetaFromEta ( 3.5 ),ThetaFromEta ( 1 ));
229  Smear::Device HcalFwd(Smear::kE, "sqrt(pow( 0.0*E, 2) + pow ( 0.5,2) *E)");
230  HcalFwd.Accept.AddZone(HcalFwdZone);
232  det.AddDevice(HcalFwd);
233 
234  // Not covered:
235  // Tracker material budget X/X0 <~5%
236  // Low-Q^2 tagger: -6.9<eta<-5.8: Delta_theta/theta < 1.5%; 10^-6 < Q2 < 10^-2 GeV2
237  // Proton spectrometer: eta>6.2: sigma_intrinsic(|t|)/|t| < 1%; Acceptance: 0.2 < pT < 1.2 GeV/c
238  // Barrel vertexing: sigma_xyz ~ 20 microns, d0(z) ~ d0(r phi) ~ (20 microns)/(pT [GeV]) + 5 microns
239  // Central muon detection
240 
241  Smear::TofBarrelSmearer tofBarrel(100, -1.0, 1.0, 10);
242  // Restricting to the tracker barrel defined above
243  // but additional cuts are done by the detector!
244  tofBarrel.Accept.AddZone( TrackBarrelZone );
245  tofBarrel.Accept.SetCharge(Smear::kCharged);
246  det.AddDevice(tofBarrel);
247 
248  return det;
249 }
250 
251 // -------------------------------------------------------------------
252 double ThetaFromEta( const double eta ) {
253  if ( !std::isnan(eta) && !std::isinf(eta) ) {
254  return 2.0 * atan( exp( -eta ));
255  }
256  throw std::runtime_error("ThetaFromEta called with NaN or Inf");
257  return -1;
258 }