eic-smear  1.0.3
A collection of ROOT classes for Monte Carlo events and a fast-smearing code simulating detector effects for the Electron-Ion Collider task force
eSTARDetector.cpp
Go to the documentation of this file.
1 
15 double etaToTheta(double eta) {
16  return 2. * std::atan(std::exp(-eta));
17 }
18 
22 Smear::Acceptance::Zone zoneEta(double etamin, double etamax) {
23  // First two arguments to Zone constructor are (theta-min, theta-max)
24  // Note that eta-min --> theta-max and eta-max --> theta-min
25  return Smear::Acceptance::Zone(etaToTheta(etamax), etaToTheta(etamin));
26 }
27 
39  // Recall that the parameterisations passed to Smear::Device objects
40  // give sigma(X) *not* sigma(X) / X
41  // Electromagnetic calorimeter in the electron (east) direction, -4 < eta < -2
42  // From Zhangbu's talk: sigma(E) / E = 2% / sqrt(E) \oplus 0.75%
43  Smear::Device eastEcal("E", "0.02 * sqrt(E) + 0.0075 * E",
45  eastEcal.Accept.AddZone(zoneEta(-4., -2.));
46  // Tracking in the electron (east) direction, -2 < eta < -1
47  // From Zhangbu's talk: sigma(p)/p = 1/(pT/pZ-1/6) x (0.45%pT \oplus 0.3%)
48  // \oplus (pZ/pT)x(0.2%/p/beta)
49  // beta = p/E so 0.2%/p/beta = 0.2%E/p^2
50  Smear::Device eastTracking("P",
51  "P/(pT/pZ-1/6)*(0.0045*pT+0.003)+pZ/pT*0.002*E/P");
52  eastTracking.Accept.AddZone(zoneEta(-2., -1.));
53  eastTracking.Accept.SetCharge(Smear::kCharged);
54  // Electromagnetic calorimeter in the barrel region, -1 < eta < 1
55  // From Zhangbu's talk: sigma(E)/E = 14%/sqrt(E) \oplus 2%
56  Smear::Device barrelEcal("E", "0.14*sqrt(E)+0.02*E", Smear::kElectromagnetic);
57  barrelEcal.Accept.AddZone(zoneEta(-1., 1.));
58  // Tracking in the barrel region, -1 < eta < 1
59  // From Zhangbu's talk: sigma(p)/p = 0.45%pT \oplus 0.3% \oplus 0.2%/p/beta
60  Smear::Device barrelTracking("P", "0.0045*pT*P+0.003*P+0.002*E/P");
61  barrelTracking.Accept.AddZone(zoneEta(-1., 1.));
62  barrelTracking.Accept.SetCharge(Smear::kCharged);
63  // Tracking in hadron (west) direction, 1 < eta < 1.7
64  // From Zhangbu's talk: as for east tracking, except 1/4 instead of 1/6
65  Smear::Device westTracking("P",
66  "P/(pT/pZ-1/4)*(0.0045*pT+0.003)+pZ/pT*0.002*E/P");
67  westTracking.Accept.SetCharge(Smear::kCharged);
68  westTracking.Accept.AddZone(zoneEta(1., 1.7));
69  // Electromagnetic calorimeter in the hadron (west) direction, 1 < eta < 2
70  // From Zhangbu's talk: sigma(E)/E = 16%/sqrt(E) \oplus 2%
71  Smear::Device westEcal("E", "0.16*sqrt(E)+0.02*E", Smear::kElectromagnetic);
72  westEcal.Accept.AddZone(zoneEta(1., 2.));
73  // Electromagnetic calorimeter in the far-hadron (west) direction
74  // 2.5 < eta < 5
75  // From Zhangbu's talk: sigma(E)/E = 12%/sqrt(E) \oplus 1.4%
76  Smear::Device westEcal2("E", "0.12*sqrt(E)+0.014*E", Smear::kElectromagnetic);
77  westEcal2.Accept.AddZone(zoneEta(2.5, 5.));
78  // Hadronic calorimeter in the far-hadron (west) direction, 2.5 < eta < 5
79  // From Zhangbu's talk: sigma(E)/E = 38%/sqrt(E) \oplus 3%
80  Smear::Device westHcal("E", "0.38*sqrt(E)+0.03*E", Smear::kHadronic);
81  westHcal.Accept.AddZone(zoneEta(2.5, 5));
82  // Assume perfect theta and phi performance i.e. momentum resolution
83  // dominates over angular resolution
84  Smear::Device theta("theta", "0");
85  Smear::Device phi("phi", "0");
86  // PID performance is unparameterised as of now
87  Smear::PerfectID pid;
88  // Combine the devices into a detector.
89  Smear::Detector estar;
90  estar.AddDevice(eastEcal);
91  estar.AddDevice(barrelEcal);
92  estar.AddDevice(westEcal);
93  estar.AddDevice(westEcal2);
94  estar.AddDevice(eastTracking);
95  estar.AddDevice(barrelTracking);
96  estar.AddDevice(westTracking);
97  estar.AddDevice(westHcal);
98  estar.AddDevice(theta);
99  estar.AddDevice(phi);
100  estar.AddDevice(pid);
101  estar.SetEventKinematicsCalculator("NM JB DA");
102  return estar;
103 }
etaToTheta
double etaToTheta(double eta)
Helper function to convert eta to theta (radians)
Definition: eSTARDetector.cpp:15
Smear::Smearer::Accept
Acceptance Accept
Definition: Smearer.h:50
Smear::Acceptance::Zone
A single contiguous region of acceptance.
Definition: Acceptance.h:62
Smear::Acceptance::SetCharge
void SetCharge(ECharge charge)
Select the charges of particles to accept.
Definition: Acceptance.cxx:38
Smear::Detector
The detector structure.
Definition: Detector.h:44
Smear::kCharged
@ kCharged
Definition: Smear.h:54
zoneEta
Smear::Acceptance::Zone zoneEta(double etamin, double etamax)
Helper function producing a zone in eta.
Definition: eSTARDetector.cpp:22
BuildDetector
Smear::Detector BuildDetector()
Smearing parameterisations for the eSTAR detector.
Definition: eSTARDetector.cpp:38
Smear::Device
Performs smearing of a single kinematic variable according to a simple expression defined via a strin...
Definition: Device.h:44
Smear::Detector::SetEventKinematicsCalculator
void SetEventKinematicsCalculator(TString)
Set the method for calculating event kinematics if FillEventKinematics is used.
Definition: Detector.cxx:67
Smear::Detector::AddDevice
void AddDevice(Smearer &device)
Adds a copy of the smearing device to this detector.
Definition: Detector.cxx:63
Smear::PerfectID
Smearer that copies the PDG ID of a particle to a smeared particle with no modification.
Definition: PerfectID.h:27
Smear::Acceptance::AddZone
void AddZone(const Zone &)
Add a new zone with user-specified coverage.
Definition: Acceptance.cxx:26
Smear::kHadronic
@ kHadronic
Definition: Smear.h:49
Smear::kElectromagnetic
@ kElectromagnetic
Definition: Smear.h:49