KaliVeda
Toolkit for HIC analysis
KVElasticScatter.cpp
1 /*
2 $Id: KVElasticScatter.cpp,v 1.9 2007/04/04 10:39:17 ebonnet Exp $
3 $Revision: 1.9 $
4 $Date: 2007/04/04 10:39:17 $
5 */
6 
7 //Created by KVClassFactory on Thu Oct 19 22:31:02 2006
8 //Author: John Frankland
9 
10 #include "KVElasticScatter.h"
11 #include "KVMultiDetArray.h"
12 #include "TH1F.h"
13 #include "KVDetector.h"
14 #include "KVGroup.h"
15 #include "KVTarget.h"
16 #include "KV2Body.h"
17 #include "TObjArray.h"
18 #include "KVExpDB.h"
19 #include "TRandom.h"
20 using namespace std;
21 
23 
24 
25 
29 {
30  //Default constructor
31  if (!gMultiDetArray)
32  Error("KVElasticScatter", "First define the multidetector array for the simulation");
33  SetRun(run);
34  gMultiDetArray->SetSimMode(kTRUE);
35 }
36 
37 
38 
39 
42 
43 KVElasticScatter::~KVElasticScatter()
44 {
45  //Destructor
46  if (fDepth)
47  delete fDepth;
48  if (fTheta)
49  delete fTheta;
50  gMultiDetArray->SetSimMode(kFALSE);
51 }
52 
53 
54 
55 
62 
64 {
65  //Set beam, detector parameters, target, etc. for run
66  //
67  // The default reaction kinematics is that of the KVDBSystem associated to the run.
68  //
69  // The default scattering position in the target is random.
70  gMultiDetArray->SetParameters(run);
71  fTarget = gMultiDetArray->GetTarget();
72  fTarget->SetRandomized();
73  fAtomicDensity = fTarget->GetAtomsPerCM2() * 1.e-24; //in barn^-1 units
74  fIntLayer = 0;
75  fMultiLayer = (fTarget->NumberOfLayers() > 1);
76  // get reaction propoerties from system
77  auto sys = gExpDB->GetDBRun(run)->GetSystem();
78  fKinematics = sys->GetKinematics();
79 }
80 
81 
82 
83 
90 
92 {
93  // Give name of detector towards which elastic scattering will occur.
94  //
95  // \note this does not necessarily mean that the particle will actually
96  // arrive in the given detector (it may stop before it), but the aim is to
97  // define the trajectory of the particle in the array.
98 
99  fDetector = gMultiDetArray->GetDetector(det);
100  if (!fDetector) {
101  Error("SetDetector", "Detector %s is unknown!", det);
102  return;
103  }
104  // get trajectory for detector
105  if (fDetector->GetNode()->GetNTraj() > 1) {
106  Warning("SetDetector", "Ambiguous trajectory choice: detector %s has %d trajectories", det, fDetector->GetNode()->GetNTraj());
107  }
108  auto traj = (KVGeoDNTrajectory*)fDetector->GetNode()->GetTrajectories()->First();
109  Info("SetDetector", "Using trajectory: %s", traj->GetName());
110 
111  //we store the association between detector type and index in list
112  fDetInd.Clear();
113  fAlignedDetectors.Clear();
114 
115  traj->IterateBackFrom();
116 
117  Int_t i = 0;
118  KVGeoDetectorNode* gdn = nullptr;
119  while ((gdn = traj->GetNextNode())) {
120  KVDetector* d = gdn->GetDetector();
121  fAlignedDetectors.Add(d);
122  fDetInd.SetValue(d->GetLabel(), i);
123  ++i;
124  }
125  //store number of aligned detectors
126  fNDets = i;
127  fHistos.Expand(fNDets);
128  fHistos.SetOwner(); //will delete histos it stores
129 }
130 
131 
132 
133 
140 
142 {
143  //For multilayer targets, use this method to choose in which
144  //layer the scattering will take place.
145  //
146  //If name="", reset any previous choice so that scattering
147  //can take place in any layer
148  if (!fTarget) {
149  cout <<
150  "<KVElasticScatter::SetTargetScatteringLayer> : No target set. Set run first."
151  << endl;
152  return;
153  }
154  fIntLayer = fTarget->GetLayerIndex(name);
155  if (fIntLayer)
156  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
157 }
158 
159 
160 
161 
165 
167 {
168  //Set binning of the GetEnergy histogram
169  // Default value is 500
170  fBinE = nbins;
171 }
172 
173 
174 
176 
178 {
179  if (!fDetector) {
180  cout <<
181  "<KVElasticScatter::CalculateScattering> : Set detector first" <<
182  endl;
183  return false;
184  }
185  fHistos.Clear(); //delete any previous histograms
186 
187  KVDetector* d;
188  TIter n(&fAlignedDetectors);
189  while ((d = (KVDetector*) n()))
190  fHistos.Add(new TH1F(Form("hEloss_%s", d->GetName()), "Eloss (MeV)", fBinE, 0., 0.));
191 
192  return true;
193 }
194 
195 
196 
198 
199 std::pair<double, double> KVElasticScatter::get_random_angles_for_scattering(const KV2Body& scattering_kinematics)
200 {
201  auto max_scattering_angle = scattering_kinematics.GetMaxAngleLab(fNucleusOfInterest);
202  if (fDetector->GetEntranceWindow().GetThetaMin() > max_scattering_angle) {
203  Warning("CalculateScattering",
204  "Detector %s has minimum theta %f deg. > max scattering angle of %s, %f deg.: abandon calculation",
205  fDetector->GetName(), fDetector->GetEntranceWindow().GetThetaMin(), scattering_kinematics.GetNucleus(fNucleusOfInterest)->GetSymbol(),
206  max_scattering_angle);
207  return {-1., -1.};
208  }
209  //set random direction of outgoing projectile
210  double th, ph;
211  th = ph = 0.;
212  do {
213  fDetector->GetRandomAngles(th, ph);
214  }
215  while (th > max_scattering_angle);
216 
217  return {th, ph};
218 }
219 
220 
221 
223 
224 void KVElasticScatter::detect_particle_fill_histograms(KVNucleus* ejectile, double theta, double phi, double xsec)
225 {
226  int j = 0;
227  TIter n(&fAlignedDetectors);
228  KVDetector* d;
229  while ((d = (KVDetector*) n())) {
230  auto cew = d->GetCentreOfEntranceWindow();
231  d->DetectParticle(ejectile, &cew);
232  //fill histograms
233  ((TH1F*)fHistos[j++])->Fill(d->GetEnergy(), xsec * sin(theta * TMath::DegToRad()));
234  //prepare for next round: set energy loss to zero
235  d->Clear();
236  }
237 }
238 
239 
240 
243 
245 {
246  // Perform scattering 'N' times
247 
248  if (!initial_checks_and_reset())
249  return;
250 
251  fNtirages = N;
252 
253  //set up histograms
254  // -- histograms for debugging: depth in target and angular distribution
255  if (fDepth)
256  delete fDepth;
257  if (fTheta)
258  delete fTheta;
259  fDepth =
260  new TH1F("hDepth", "Depth (mg/cm2)", 500, 0., fTarget->GetTotalEffectiveThickness());
261  fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
262 
263  //set random interaction point for scattering
264  TVector3 IP;
265  if (fIntLayer) {
266  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
267  IP = fTarget->GetInteractionPoint();
268  }
269  else {
270  IP = fTarget->GetInteractionPoint(fBeamDirection);
271  }
272 
273  /* -------------------------------------------------------------------------------------------------------------------------- */
274 
275  if (fGetTargFromLay) {
276  // get target nucleus properties from scattering layer ?
277  KVMaterial* targ_mat = fTarget->GetLayer(IP);
278  auto zt = targ_mat->GetZ();
279  auto at = targ_mat->GetMass();
280  // check for fractional Z or A - compound material!
281  if ((zt != int(zt)) || at != int(at)) {
282  Error("CalculateScattering", "Target has non-elemental layer with effective Z=%f and A=%f. Cannot use.",
283  zt, at);
284  fGetTargFromLay = false;
285  }
286  else
287  fKinematics->SetTarget({int(zt), int(at)});
288  }
289 
290  for (int i = 0; i < N; i++) {
291 
292  reset_before_new_scattering();
293 
294  auto Projectile = *(fKinematics->GetNucleus(KV2Body::projectile));
295  auto Target = *(fKinematics->GetNucleus(KV2Body::target));
296 
297  //calculate slowing of incoming projectile
298  fTarget->SetIncoming();
299  fTarget->DetectParticle(&Projectile);
300 
301  KV2Body scattering{Projectile, Target, fExx};
302  if (fOutgoingPL)
303  scattering.SetOutgoing(*fOutgoingPL);
304 
305  if (fExx > 0)
306  scattering.SetEDiss(fExx);
307 
308  scattering.CalculateKinematics();
309 
310  auto theta_phi = get_random_angles_for_scattering(scattering);
311 
312  //set energy of scattered nucleus
313  auto elab = scattering.GetELab(fNucleusOfInterest, theta_phi.first);
314  auto ejectile = scattering.GetNucleus(fNucleusOfInterest);
315 
316  auto kine_branch = fKineSol;
317  if (fKineSol == KV2Body::low_E_branch) {
318  if (!elab->second()) {
319  Warning("CalculateScattering", "Low energy branch requested, but only high energy kinematics for theta = %f deg.",
320  theta_phi.first);
321  ejectile->SetEnergy(*elab->first());
322  kine_branch = KV2Body::high_E_branch;
323  }
324  else {
325  ejectile->SetEnergy(*elab->second());
326  }
327  }
328  else if (fKineSol == KV2Body::both_branches) {
329  if (!elab->second()) {
330  Warning("CalculateScattering", "Both branches requested, but only high energy kinematics for theta = %f deg.",
331  theta_phi.first);
332  ejectile->SetEnergy(*elab->first());
333  kine_branch = KV2Body::high_E_branch;
334  }
335  else {
336  ejectile->SetEnergy(1.0);
337  ejectile->SetTheta(theta_phi.first);
338  auto angle_with_beam = TMath::RadToDeg() * ejectile->Angle(fBeamDirection);
339  auto xsec_hi = TMath::Abs(scattering.GetXSecRuthLab(angle_with_beam, fNucleusOfInterest, KV2Body::high_E_branch));
340  auto xsec_lo = TMath::Abs(scattering.GetXSecRuthLab(angle_with_beam, fNucleusOfInterest, KV2Body::low_E_branch));
341  if (gRandom->Uniform() < xsec_hi / (xsec_hi + xsec_lo)) {
342  // high
343  ejectile->SetEnergy(*elab->first());
344  kine_branch = KV2Body::high_E_branch;
345  }
346  else {
347  // low
348  ejectile->SetEnergy(*elab->second());
349  kine_branch = KV2Body::low_E_branch;
350  }
351  }
352  }
353  else
354  ejectile->SetEnergy(*elab->first());
355 
356  ejectile->SetTheta(theta_phi.first);
357  ejectile->SetPhi(theta_phi.second);
358 
359  // get Rutherford Xsec for projectile scattered to theta degrees in lab
360  // taking into account the initial direction of the beam (if beam direction = (0,0,1) this is just the same as theta)
361  double angle_with_beam = TMath::RadToDeg() * ejectile->Angle(fBeamDirection);
362  auto xsec = TMath::Abs(scattering.GetXSecRuthLab(angle_with_beam, fNucleusOfInterest, kine_branch));
363 
364  // fill angular distribution
365  fTheta->Fill(theta_phi.first, xsec);
366 
367  //slowing of outgoing ejectile in target
368  fTarget->SetOutgoing();
369 
370  fTarget->DetectParticle(ejectile);
371  fDepth->Fill(IP.z());
372 
373  //now detect particle in detector(s)
374  detect_particle_fill_histograms(ejectile, theta_phi.first, theta_phi.second, xsec);
375 
376  //set random interaction point for scattering
377  if (fIntLayer) {
378  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
379  IP = fTarget->GetInteractionPoint();
380  }
381  else {
382  IP = fTarget->GetInteractionPoint(fBeamDirection);
383  }
384  if (fGetTargFromLay) {
385  // get target nucleus properties from scattering layer ?
386  KVMaterial* targ_mat = fTarget->GetLayer(IP);
387  auto zt = targ_mat->GetZ();
388  auto at = targ_mat->GetMass();
389  // check for fractional Z or A - compound material!
390  if ((zt != int(zt)) || at != int(at)) {
391  Error("CalculateScattering", "Target has non-elemental layer with effective Z=%f and A=%f. Cannot use.",
392  zt, at);
393  fGetTargFromLay = false;
394  }
395  else
396  fKinematics->SetTarget({int(zt), int(at)});
397  }
398  }
399  end_of_run();
400 }
401 
402 
403 
404 
405 
408 
410 {
411  //Energy loss in detector with given 'label' through which scattered particle passes.
412 
413  return (fDetInd.HasParameter(label) ? GetEnergy(fDetInd.GetIntValue(label)) : 0);
414 }
415 
416 
417 
418 
425 
427 {
428  //Energy loss in any detector through which scattered particle passes.
429  //
430  //The index corresponds to the order in which detectors are crossed by the
431  //particle, beginning with 0 for the first detector, and ending with
432  //(GetNDets()-1)
433  return (TH1F*)fHistos[index];
434 }
435 
436 
437 
440 
442 {
443  // Print details of calculation to be performed
444 
445  Info("Print", "CALCULATION DETAILS:");
446  if (fGetTargFromLay)
447  Info("Print", "Projectile is %s %.1f MeV/u; target nucleus depends on scattering layer.",
448  fKinematics->GetNucleus(KV2Body::projectile)->GetSymbol(),
449  fKinematics->GetNucleus(KV2Body::projectile)->GetAMeV());
450  else
451  Info("Print", "Entrance channel is: %s", fKinematics->GetTitle());
452  if (fExx > 0) {
453  Info("Print", "Inelastic scattering: outgoing Ex=%f", fExx);
454  }
455  Info("Print", "Target used is:");
456  fTarget->Print();
457  Info("Print", "Scattering will be calculated for outgoing %s",
458  fNucleusOfInterest == KV2Body::projectile_like ? "PROJECTILE_LIKE" : "TARGET_LIKE");
459  Info("Print", "We will use %s kinematic solutions",
460  fKineSol == KV2Body::both_branches ? "BOTH" : (fKineSol == KV2Body::high_E_branch ? "the HIGH ENERGY" : "the LOW ENERGY"));
461  if (fDetector) {
462  Info("Print", "Particles scattered towards %s: theta:[%f,%f] phi:[%f,%f]", fDetector->GetName(),
463  fDetector->GetThetaMin(), fDetector->GetThetaMax(), fDetector->GetPhiMin(), fDetector->GetPhiMax());
464  }
465 }
466 
467 
468 //__________________________________________________________________//
469 
int Int_t
#define d(i)
char Char_t
constexpr Bool_t kTRUE
#define N
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t index
char name[80]
R__EXTERN TRandom * gRandom
char * Form(const char *fmt,...)
Relativistic binary kinematics calculator.
Definition: KV2Body.h:178
Double_t GetMaxAngleLab(nucleus_of_interest i) const
Definition: KV2Body.cpp:546
void SetOutgoing(const KVNucleus &proj_out)
Definition: KV2Body.cpp:411
@ target
Definition: KV2Body.h:265
@ projectile_like
Definition: KV2Body.h:266
@ projectile
Definition: KV2Body.h:264
@ low_E_branch
Definition: KV2Body.h:271
@ both_branches
Definition: KV2Body.h:272
@ high_E_branch
Definition: KV2Body.h:270
KVNucleus * GetNucleus(nucleus_of_interest i) const
Definition: KV2Body.cpp:468
KVDBSystem * GetSystem() const
Definition: KVDBRun.cpp:240
KV2Body * GetKinematics()
Definition: KVDBSystem.cpp:103
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:160
Calculate elastic scattering spectra in specific detectors of a multidetector array ,...
virtual std::pair< double, double > get_random_angles_for_scattering(const KV2Body &scattering_kinematics)
void SetDetector(const Char_t *det)
virtual bool initial_checks_and_reset()
TH1F * GetEnergy()
Return pointer to energy loss histogram for chosen detector (in MeV)
void CalculateScattering(Int_t N)
Perform scattering 'N' times.
void Print() const
Print details of calculation to be performed.
void SetRun(Int_t run)
void SetTargetScatteringLayer(const Char_t *name)
virtual void detect_particle_fill_histograms(KVNucleus *ejectile, double theta, double phi, double xsec)
void SetEbinning(Int_t nbins)
KVDBRun * GetDBRun(Int_t number) const
Definition: KVExpDB.h:144
Path taken by particles through multidetector geometry.
Information on relative positions of detectors & particle trajectories.
KVDetector * GetDetector() const
virtual KVDetector * GetDetector(const Char_t *name) const
Return detector in this structure with given name.
Description of physical materials used to construct detectors & targets; interface to range tables.
Definition: KVMaterial.h:90
Double_t GetZ() const
Definition: KVMaterial.cpp:390
Double_t GetMass() const
Definition: KVMaterial.cpp:302
KVTarget * GetTarget()
virtual void SetParameters(UInt_t n, Bool_t physics_parameters_only=kFALSE)
virtual void SetSimMode(Bool_t on=kTRUE)
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
const Char_t * GetSymbol(Option_t *opt="") const
Definition: KVNucleus.cpp:71
void SetRandomized(Bool_t r=kTRUE)
Definition: KVTarget.h:254
virtual Double_t Uniform(Double_t x1, Double_t x2)
Double_t z() const
RVec< PromoteType< T > > sin(const RVec< T > &v)
const Int_t n
void Error(const char *location, const char *fmt,...)
void Info(const char *location, const char *fmt,...)
void Warning(const char *location, const char *fmt,...)
constexpr Double_t DegToRad()
Double_t Abs(Double_t d)
constexpr Double_t RadToDeg()
ClassImp(TPyArg)