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
ePHENIXDetector.cpp
Go to the documentation of this file.
1 
10 /*
11  Example smearing script for the ePHENIX detector.
12 
13  It defines:
14  - EPhenixMomentum: a special smearing class for the ePHENIX momentum
15  performance
16  - BuildDetector: a function to generate the full ePHENIX detector description
17 
18  This script must be compiled in ROOT before running.
19  Therefore, you must first make sure that
20  (1) libeicsmear is loaded, via
21  gSystem->Load("libeicsmear");
22  (2) the eicsmear headers are accessible, by doing e.g.
23  gSystem->AddIncludePath(" -I/path/to/eic-smear/include");
24  You can then compile it in a ROOT session via
25  .L ePHENIXDetector.cpp+g
26  */
27 
28 #include <algorithm> // For std::max
29 #include <cmath>
30 #include <list>
31 
32 #include <TCollection.h> // For TIter
33 #include <TGraph.h>
34 #include <TList.h>
35 #include <TLorentzVector.h>
36 #include <TMultiGraph.h>
37 #include <TPad.h>
38 #include <TRandom.h>
39 
42 #include "eicsmear/smear/Device.h"
44 #include "eicsmear/smear/Smearer.h"
47 
58  public:
62  virtual ~EPhenixMomentum();
70  EPhenixMomentum(bool multipleScattering = true);
77  void Initialise();
83  TMultiGraph* Graphs();
87  virtual EPhenixMomentum* Clone(const char* /* unused */) const;
91  virtual void Smear(const erhic::VirtualParticle& particle,
92  Smear::ParticleMCS& smeared);
98  virtual double computeMultipleScattering(
99  const erhic::VirtualParticle& particle) const;
103  virtual void Draw(Option_t* option = "ac");
104 
105  private:
106  TMultiGraph* mGraphDrawer;
107  Bool_t mMultipleScattering;
108  ClassDef(EPhenixMomentum, 1)
109 };
110 
111 EPhenixMomentum::EPhenixMomentum(bool multipleScattering)
112  : mGraphDrawer(NULL), mMultipleScattering(multipleScattering) {
113  // Only use charged particles, as this is a tracker
115 }
116 
118  if (mGraphDrawer) {
119  delete mGraphDrawer;
120  mGraphDrawer = NULL;
121  } // if
122 }
123 
125  // Pseudorapidity values
126  const double eta[75] = {
127  -3.0, -2.9, -2.8, -2.7, -2.6, -2.5, -2.4, -2.3, -2.2, -2.1, -2.0, -2.0,
128  -1.9, -1.8, -1.7, -1.6, -1.5, -1.5, -1.4, -1.3, -1.2, -1.1, -1.0, -0.9,
129  -0.8, -0.7, -0.6, -0.5, -0.4, -0.3, -0.2, -0.1, 0.0, 0.1, 0.2, 0.3,
130  0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1.0, 1.0, 1.1, 1.2, 1.3, 1.4,
131  1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.1, 2.2, 2.3, 2.4, 2.5, 2.5,
132  2.6, 2.7, 2.8, 2.9, 3.0, 3.1, 3.2, 3.3, 3.4, 3.5, 3.6, 3.7,
133  3.8, 3.9, 4.0
134  };
135  // sigma(1/p) values at each value of eta
136  const double sigma[75] = {
137  0.010500, 0.009820, 0.009140, 0.008460, 0.007780, 0.007100, 0.006520,
138  0.005940, 0.005360, 0.004780, 0.004200, 0.009811, 0.008759, 0.007479,
139  0.006242, 0.005436, 0.004524, 0.006887, 0.005374, 0.004485, 0.003822,
140  0.003164, 0.002592, 0.002791, 0.002991, 0.003187, 0.003374, 0.003547,
141  0.003700, 0.003827, 0.003921, 0.003980, 0.004000, 0.003980, 0.003921,
142  0.003827, 0.003700, 0.003547, 0.003374, 0.003187, 0.002991, 0.002791,
143  0.002592, 0.001200, 0.001280, 0.001360, 0.001440, 0.001520, 0.001600,
144  0.001840, 0.002080, 0.002320, 0.002560, 0.002800, 0.003160, 0.003520,
145  0.003880, 0.004240, 0.004600, 0.002300, 0.002620, 0.002940, 0.003260,
146  0.003580, 0.003900, 0.004400, 0.004900, 0.005400, 0.005900, 0.006400,
147  0.007220, 0.008040, 0.008860, 0.009680, 0.010500
148  };
149  mGraphDrawer = new TMultiGraph;
150  mGraphDrawer->SetTitle(";#eta;Momentum Resolution, #delta(1/p) (1/GeV)");
151  // Use different graphs for each region to avoid problems at discontinuities
152  // both with the calculated values and with drawing
153  int points[5] = {11, 6, 26, 16, 16}; // Number of points in each graph
154  int colours[5] = {kBlue, kMagenta, kRed, kGreen + 3, kGreen - 2};
155  int position(0); // Tracks cumulative array index
156  for (int i(0); i < 5; ++i) {
157  TGraph* graph = new TGraph(points[i], eta + position, sigma + position);
158  TString name("graph");
159  name += i; // Append graph number to name
160  graph->SetName(name);
161  // Apply formatting
162  graph->SetLineWidth(3);
163  graph->SetLineColor(colours[i]);
164  graph->SetMarkerColor(colours[i]);
165  mGraphDrawer->Add(graph);
166  position += points[i]; // Update total array index
167  } // for
168 }
169 
170 TMultiGraph* EPhenixMomentum::Graphs() {
171  if (!mGraphDrawer) {
172  Initialise();
173  } // if
174  return mGraphDrawer;
175 }
176 
178  EPhenixMomentum* clone = new EPhenixMomentum(*this);
179  // We need to "de-initialise" the clone so it doesn't retain the
180  // same graphs as this instance.
181  clone->mGraphDrawer = NULL;
182  return clone;
183 }
184 
186  Smear::ParticleMCS& smeared) {
187  TGraph* graph(NULL);
188  TGraph* iterGraph(NULL); // For iterator position
189  TIter iter(Graphs()->GetListOfGraphs());
190  // Loop through graphs for each eta range until we find the one for
191  // the eta of this particle
192  while ((iterGraph = static_cast<TGraph*>(iter()))) {
193  // Array of eta points in the graph
194  double* eta = iterGraph->GetX();
195  if (particle.GetEta() > eta[0] &&
196  particle.GetEta() < eta[iterGraph->GetN() - 1] ) {
197  graph = iterGraph;
198  break;
199  } // if
200  } // while
201  if (graph) {
202  // Stored values are sigma(1/p).
203  // We want to get sigma(p), which is equivalent to sigma(1/p) * p^2.
204  double sigmaP = graph->Eval(particle.GetEta()) // sigma(1/p)
205  * std::pow(particle.GetP(), 2.);
206  // Add multiple scattering in quadrature with the linear term if requested
207  if (mMultipleScattering) {
208  sigmaP = std::sqrt(std::pow(sigmaP, 2.) +
209  std::pow(computeMultipleScattering(particle), 2.));
210  } // if
211  double sigma = gRandom->Gaus(particle.GetP(), sigmaP);
212  smeared.SetP(std::max(sigma, 0.)); // Store p, must not be negative
213  } // if
214 }
215 
217  const erhic::VirtualParticle& particle) const {
218  // Here is a duplicate of the relevant information from Sasha's email:
219  // "As for multiple scattering term, what has been simulated so far is eta
220  // ranges 2-3 and 3-4. Numbers keep changing after inclusion of different
221  // material (associated with detectors and read-out). The current very
222  // conservative advice is to use 3% for eta=2-3; in eta=3-4, log of this
223  // term changes ~linearly from 3% at eta=3 to ~15% at eta=4."
224  // Note that the 3% etc is sigma(P) / P.
225  const double eta = particle.GetEta();
226  if (eta > 2. && eta <= 3.) {
227  return 0.03 * particle.GetP(); // Flat sigma(P)/P = 3% in this range
228  } else if (eta > 3. && eta <= 4.) {
229  // The behaviour Sasha describes is (I think!)
230  // log(M.S.) ~= log(3%) + (log(15%) - log(3%)) * (eta - 3) [3 < eta < 4]
231  double logSigmaP = std::log(0.03) +
232  (std::log(0.15) - std::log(0.03)) * (eta - 3);
233  // Remember to multiply by momentum, as the 3-15% is sigma(P)/P and we
234  // want sigma(P)
235  return std::exp(logSigmaP) * particle.GetP();
236  } // if
237  return 0.; // Outside 2 < eta < 4
238 }
239 
240 void EPhenixMomentum::Draw(Option_t* option) {
241  Graphs()->Draw(option);
242  if (gPad) {
243  gPad->SetLogy(1);
244  } // if
245 }
246 
252 double etaToTheta(double eta) {
253  return 2. * std::atan(std::exp(-eta));
254 }
255 
270 Smear::Detector BuildDetector(bool multipleScattering = true) {
271  EPhenixMomentum momentum(multipleScattering);
272  // Define acceptance zones for different ePHENIX regions:
273  // - electron-going: -4 < eta < -1
274  // - barrel: -1 < eta < 1
275  // - hadron-going: 1 < eta < 4
276  // As the same calorimeter performance is used in the barrel and hadron-going
277  // directions we define them as a single zone
278  // Note etamin -> thetamax, etamax -> thetamin
279  Smear::Acceptance::Zone electronDirection(etaToTheta(-1.), etaToTheta(-4.));
280  Smear::Acceptance::Zone barrelAndHadronDirection(etaToTheta(4.),
281  etaToTheta(-1.));
282  // Electron-going electromagnetic calorimeter
283  Smear::Device electronEcal("E", "0.015*sqrt(E) + 0.01*E",
285  electronEcal.Accept.AddZone(electronDirection);
286  // Barrel and hadron-going electromagnetic calorimeter
287  Smear::Device barrelAndHadronEcal("E", "0.12*sqrt(E) + 0.02*E",
289  barrelAndHadronEcal.Accept.AddZone(barrelAndHadronDirection);
290  // Barrel and hadron-going hadronic calorimeter
291  Smear::Device barrelAndHadronHcal("E", "sqrt(E)", Smear::kHadronic);
292  barrelAndHadronHcal.Accept.AddZone(barrelAndHadronDirection);
293  // Assume perfect theta and phi performance i.e. momentum resolution
294  // dominates over angular resolution
295  Smear::Device theta("theta", "0");
296  Smear::Device phi("phi", "0");
297  // PID performance is unparameterised as of now
298  Smear::PerfectID pid;
299  // Combine the devices into a detector.
300  Smear::Detector ephenix;
301  ephenix.AddDevice(momentum);
302  ephenix.AddDevice(electronEcal);
303  ephenix.AddDevice(barrelAndHadronEcal);
304  ephenix.AddDevice(barrelAndHadronHcal);
305  ephenix.AddDevice(theta);
306  ephenix.AddDevice(phi);
307  ephenix.AddDevice(pid);
308  ephenix.SetEventKinematicsCalculator("NM JB DA");
309  return ephenix;
310 }
311 
312 /*
313 Here is the email from Sasha giving the calorimeter performances used above:
314 
315 From: Alexander Bazilevsky <shura@bnl.gov>
316  Subject: Re: a favour
317  Date: 26 December, 2013 1:47:19 PM EST
318  To: elke-caroline aschenauer <elke@bnl.gov>, Thomas P Burton <tpb@bnl.gov>
319 
320 Hello Elke, Thomas,
321 
322 Of course numbers are still drifting, but for now what you could use is roughly the following:
323 
324 EMCal sigma_E/E:
325 e-going direction: 1.5%/sqrt(E) \oplus 1%
326 h-going direction and barrel: 12%/sqrt(E) \oplus 2%
327 
328 HCal sigma_E/E:
329 h-going direction and barrel: 100%/sqrt(E)
330 
331 Tracking sigma_p/p:
332 is quite complicated: the linear term is in attached plot; we also calculated constant (multiple scattering) term from Geant, which appeared to vary from under 1% at eta~1 to 3% at eta~3 and to larger values (~10%) towards eta~4. I'll give you numerical parametrization one of the following days, as soon as I get them from the author of these calculations (Jin Huang).
333 
334 Regards,
335 
336 Sasha.
337 
338 --------------------------------------------------------------------------------
339 
340 Here is the email from Sasha giving the momentum resolution points used above:
341 
342 From: Alexander Bazilevsky <shura@bnl.gov>
343  Subject: Re: a favour
344  Date: 2 January, 2014 11:01:39 AM EST
345  To: elke-caroline aschenauer <elke@bnl.gov>, Thomas P Burton <tpb@bnl.gov>
346 
347 Hello Elke, Thomas,
348 
349 Ok, the linear term in momentum resolution (as in the plot I sent you before):
350 
351 eta d(1/p) in (1/GeV)
352 -3.0 0.010500
353 -2.9 0.009820
354 -2.8 0.009140
355 -2.7 0.008460
356 -2.6 0.007780
357 -2.5 0.007100
358 -2.4 0.006520
359 -2.3 0.005940
360 -2.2 0.005360
361 -2.1 0.004780
362 -2.0 0.004200
363 
364 -2.0 0.009811
365 -1.9 0.008759
366 -1.8 0.007479
367 -1.7 0.006242
368 -1.6 0.005436
369 -1.5 0.004524
370 
371 -1.5 0.006887
372 -1.4 0.005374
373 -1.3 0.004485
374 -1.2 0.003822
375 -1.1 0.003164
376 -1.0 0.002592
377 -0.9 0.002791
378 -0.8 0.002991
379 -0.7 0.003187
380 -0.6 0.003374
381 -0.5 0.003547
382 -0.4 0.003700
383 -0.3 0.003827
384 -0.2 0.003921
385 -0.1 0.003980
386 0.0 0.004000
387 0.1 0.003980
388 0.2 0.003921
389 0.3 0.003827
390 0.4 0.003700
391 0.5 0.003547
392 0.6 0.003374
393 0.7 0.003187
394 0.8 0.002991
395 0.9 0.002791
396 1.0 0.002592
397 
398 1.0 0.001200
399 1.1 0.001280
400 1.2 0.001360
401 1.3 0.001440
402 1.4 0.001520
403 1.5 0.001600
404 1.6 0.001840
405 1.7 0.002080
406 1.8 0.002320
407 1.9 0.002560
408 2.0 0.002800
409 2.1 0.003160
410 2.2 0.003520
411 2.3 0.003880
412 2.4 0.004240
413 2.5 0.004600
414 
415 2.5 0.002300
416 2.6 0.002620
417 2.7 0.002940
418 2.8 0.003260
419 2.9 0.003580
420 3.0 0.003900
421 3.1 0.004400
422 3.2 0.004900
423 3.3 0.005400
424 3.4 0.005900
425 3.5 0.006400
426 3.6 0.007220
427 3.7 0.008040
428 3.8 0.008860
429 3.9 0.009680
430 4.0 0.010500
431 
432 As for multiple scattering term, what has been simulated so far is eta ranges 2-3 and 3-4. Numbers keep changing after inclusion of different material (associated with detectors and read-out). The current very conservative advice is to use 3% for eta=2-3; in eta=3-4, log of this term changes ~linearly from 3% at eta=3 to ~15% at eta=4.
433 
434 That's what we have for now.
435 
436 Regards,
437 
438 Sasha.
439  */
PerfectID.h
ParticleMCS.h
erhic::VirtualParticle::GetP
virtual Double_t GetP() const =0
Returns the magnitude of 3-momentum (GeV).
Smear::Smearer::Accept
Acceptance Accept
Definition: Smearer.h:50
EPhenixMomentum::Graphs
TMultiGraph * Graphs()
Returns a pointer to the graphs.
Definition: ePHENIXDetector.cpp:170
EPhenixMomentum::EPhenixMomentum
EPhenixMomentum(bool multipleScattering=true)
Constructor.
Definition: ePHENIXDetector.cpp:111
EPhenixMomentum::Smear
virtual void Smear(const erhic::VirtualParticle &particle, Smear::ParticleMCS &smeared)
Smears the input ParticleMC and stores the results in the ParticleMCS.
Definition: ePHENIXDetector.cpp:185
Smear::ParticleMCS::SetP
virtual void SetP(Double_t)
Definition: ParticleMCS.h:233
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
EPhenixMomentum::Draw
virtual void Draw(Option_t *option="ac")
Draw all graphs.
Definition: ePHENIXDetector.cpp:240
Smear::Detector
The detector structure.
Definition: Detector.h:44
erhic::VirtualParticle::GetEta
virtual Double_t GetEta() const =0
Returns the pseudorapidity.
EPhenixMomentum::Initialise
void Initialise()
Initialise the object.
Definition: ePHENIXDetector.cpp:124
VirtualParticle.h
erhic::VirtualParticle
Abstract base class for a general particle.
Definition: VirtualParticle.h:23
Smear::kCharged
@ kCharged
Definition: Smear.h:54
Smearer.h
Device.h
Smear::Smearer
Abstract base class for objects performing smearing.
Definition: Smearer.h:33
Acceptance.h
Smear::ParticleMCS
A smeared Monte Carlo particle.
Definition: ParticleMCS.h:27
EPhenixMomentum::Clone
virtual EPhenixMomentum * Clone(const char *) const
Duplicate this object.
Definition: ePHENIXDetector.cpp:177
BuildDetector
Smear::Detector BuildDetector(bool multipleScattering=true)
Smearing parameterisations for the ePHENIX detector.
Definition: ePHENIXDetector.cpp:270
Smear::Device
Performs smearing of a single kinematic variable according to a simple expression defined via a strin...
Definition: Device.h:44
EPhenixMomentum
Smearing class describing ePHENIX momentum resolution.
Definition: ePHENIXDetector.cpp:57
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
etaToTheta
double etaToTheta(double eta)
Helper function to convert eta to theta (radians)
Definition: ePHENIXDetector.cpp:252
Smear::PerfectID
Smearer that copies the PDG ID of a particle to a smeared particle with no modification.
Definition: PerfectID.h:27
EPhenixMomentum::computeMultipleScattering
virtual double computeMultipleScattering(const erhic::VirtualParticle &particle) const
Return sigma(P) due to multiple scattering.
Definition: ePHENIXDetector.cpp:216
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
Detector.h
Smear::kElectromagnetic
@ kElectromagnetic
Definition: Smear.h:49
EPhenixMomentum::~EPhenixMomentum
virtual ~EPhenixMomentum()
Destructor.
Definition: ePHENIXDetector.cpp:117