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
STARDetector.cpp
Go to the documentation of this file.
1 
14 double etaToTheta(const double eta) {
15  TLorentzVector v;
16  v.SetPtEtaPhiM(1., eta, 0., 0.);
17  return v.Theta();
18 }
19 
23 Smear::Acceptance::Zone makeZone(double etaMin, double etaMax) {
24  // Note that we need to flip the order of arguments because the numerically
25  // larger eta corresponds to the minimum angle and the smaller eta to the
26  // maximum angle.
27  // e.g. for eta -1 to 1, -1 --> theta=2.4 and 1 --> theta=0.7
28  return Smear::Acceptance::Zone(etaToTheta(etaMax), etaToTheta(etaMin));
29 }
30 
42  // Electromagnetic calorimeter.
43  // Acceptance is valid for the BEMC (-1 < eta < 1) and EEMC (1 < eta < 2)
44  // so define a single acceptance zone -1 < eta < 2.
45  // Acceptance is defined in terms of angle (theta) not pseudorapidity (eta)
46  // so we need to convert eta -> theta.
48  "0.015*E+0.14*sqrt(E)",
50  emCal.Accept.AddZone(makeZone(-1., 2.));
51  // Tracking, comprising momentum and theta.
52  // These just span the TPC, -1 < eta < 1.
53  Smear::Device momentum(Smear::kP, "0.005*P+0.004*P*P");
55  momentum.Accept.AddZone(makeZone(-1., 1.));
56  // See comments at the end of the file for Zhangbu's email about
57  // theta resolution.
59  "sqrt(pow(3e-4, 2) + pow(9e-4 / P, 2) / sin(theta))");
61  theta.Accept.AddZone(makeZone(-1., 1.));
62  // We don't have a parameterisation for phi, so just use a
63  // device that gives perfect resolution.
64  Smear::Device phi(Smear::kPhi, "0");
66  // Combine the devices into a detector.
67  Smear::Detector star;
68  star.AddDevice(emCal);
69  star.AddDevice(momentum);
70  star.AddDevice(theta);
71  star.AddDevice(phi);
72  star.SetEventKinematicsCalculator("NM JB DA");
73  return star;
74 }
75 
76 /*
77 Here is the email from Zhangbu giving the parameterisations
78 used in the above.
79 
80 
81 Date: Thu, 2 Jun 2011 17:35:45 -0400
82 From: "Xu, Zhangbu" <xzb@bnl.gov>
83 To: "Aschenauer, Elke" <elke@bnl.gov>
84 Subject: theta angular resolution
85 
86 Hi, Elke:
87 
88 Talk to a couple of people doing the TPC calibration,
89 from drift velocity and TPC cluster uncertainties, the
90 d (theta) ~=3xe-4
91 
92 I am working on the MC simulation to get the resolution
93 due to multiple scattering, this will be momentum dependent.
94 A hand calculation of the theta angular spread due to multiple
95 scattering for STAR TPC material (mainly beam pipe and inner field
96 cage ~=0.5%X0) is
97 theta0 = 13.6MeV/beta/p*sqrt(x/X0)~=0.9MeV/p,
98 So for a 2GeV/c electron, the theta angular resolution will be
99 5e-4, close to the detector resolution.
100 
101 So it is probably reasonable to assume the theta angular resolution
102 to be: d(theta) = sqrt((3xe-4)**2+(0.9MeV/p)**2/sin(theta))
103 
104 
105 Thanks!
106 
107 Zhangbu
108 */
etaToTheta
double etaToTheta(const double eta)
Convert pseudorapidity (eta) to polar angle (theta) in radians.
Definition: STARDetector.cpp:14
Smear::kE
@ kE
Definition: Smear.h:44
Smear::Smearer::Accept
Acceptance Accept
Definition: Smearer.h:50
Smear::kTheta
@ kTheta
Definition: Smear.h:44
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
Smear::kPhi
@ kPhi
Definition: Smear.h:44
BuildDetector
Smear::Detector BuildDetector()
Smearing parameterisations for the STAR detector.
Definition: STARDetector.cpp:41
Smear::kP
@ kP
Definition: Smear.h:44
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::Acceptance::AddZone
void AddZone(const Zone &)
Add a new zone with user-specified coverage.
Definition: Acceptance.cxx:26
makeZone
Smear::Acceptance::Zone makeZone(double etaMin, double etaMax)
Returns an Acceptance::Zone spanning a range in eta.
Definition: STARDetector.cpp:23
Smear::kElectromagnetic
@ kElectromagnetic
Definition: Smear.h:49