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 
19 using namespace std;
20 
22 
23 
24 
27 KVElasticScatter::KVElasticScatter(): fBeamDirection(0, 0, 1)
28 {
29  //Default constructor
30  fDepth = 0;
31  fTheta = 0;
32  fBinE = 500;
33  fEnergy = 0;
34  fKinematics = 0;
35  fTarget = 0;
36  fIntLayer = fNDets = 0;
37  fDetector = 0;
38  fMultiLayer = kFALSE;
39  fExx = 0.;
40  fHistos = 0;
41  has_outgoing = false;
42  if (gMultiDetArray) {
43  gMultiDetArray->SetSimMode(kTRUE);
44  }
45  else {
46  Warning("KVElasticScatter", "gMultiDetArray does not refer to a valid multidetector array");
47  printf("Define it before using this class, and put it in simulation mode : gMultiDetArray->SetSimMode(kTRUE)");
48  }
49 }
50 
51 
52 
53 
56 
57 KVElasticScatter::~KVElasticScatter()
58 {
59  //Destructor
60  if (fDepth)
61  delete fDepth;
62  if (fTheta)
63  delete fTheta;
64  if (fKinematics)
65  delete fKinematics;
66  if (fHistos)
67  delete fHistos;
68  gMultiDetArray->SetSimMode(kFALSE);
69 }
70 
71 
72 
73 
76 
78 {
79  //Set detector parameters, target, etc. for run
80  gMultiDetArray->SetParameters(run);
81  fTarget = gMultiDetArray->GetTarget();
82  fTarget->SetRandomized();
83  fIntLayer = 0;
84  fMultiLayer = (fTarget->NumberOfLayers() > 1);
85 }
86 
87 
88 
89 
92 
94 {
95  //Set projectile Z and A
96 
97  fProj.SetZandA(Z, A);
98 }
99 
100 
101 
102 
105 
107 {
108  //Set energy of projectile in MeV
109 
110  fProj.SetEnergy(e);
111  fEnergy = e;
112 }
113 
114 
115 
116 
123 
125 {
126  // Give name of detector towards which elastic scattering will occur.
127  //
128  // \note this does not necessarily mean that the particle will actually
129  // arrive in the given detector (it may stop before it), but the aim is to
130  // define the trajectory of the particle in the array.
131 
132  fDetector = gMultiDetArray->GetDetector(det);
133  if (!fDetector) {
134  Error("SetDetector", "Detector %s is unknown!", det);
135  return;
136  }
137  // get trajectory for detector
138  if (fDetector->GetNode()->GetNTraj() > 1) {
139  Warning("SetDetector", "Ambiguous trajectory choice: detector %s has %d trajectories", det, fDetector->GetNode()->GetNTraj());
140  }
141  auto traj = (KVGeoDNTrajectory*)fDetector->GetNode()->GetTrajectories()->First();
142  Info("SetDetector", "Using trajectory: %s", traj->GetName());
143 
144  //we store the association between detector type and index in list
145  fDetInd.Clear();
146  fAlignedDetectors.Clear();
147 
148  traj->IterateBackFrom();
149 
150  Int_t i = 0;
151  KVGeoDetectorNode* gdn = nullptr;
152  while ((gdn = traj->GetNextNode())) {
153  KVDetector* d = gdn->GetDetector();
154  fAlignedDetectors.Add(d);
155  fDetInd.SetValue(d->GetLabel(), i);
156  ++i;
157  }
158  //store number of aligned detectors
159  fNDets = i;
160 }
161 
162 
163 
164 
170 
172 {
173  //For multilayer targets, use this method to choose in which
174  //layer the scattering will take place.
175  //If name="", reset any previous choice so that scattering
176  //can take place in any layer
177  if (!fTarget) {
178  cout <<
179  "<KVElasticScatter::SetTargetScatteringLayer> : No target set. Set run first."
180  << endl;
181  return;
182  }
183  fIntLayer = fTarget->GetLayerIndex(name);
184  if (fIntLayer)
185  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
186 }
187 
188 
189 
190 
194 
196 {
197  //Set binning of the GetEnergy histogram
198  // Default value is 500
199  fBinE = nbins;
200 }
201 
202 
203 
207 
209 {
210  //Perform scattering 'N' times for current values
211  //of particle Z, A and energy, target excited state, and detector.
212 
213  if (!fProj.IsDefined()) {
214  cout <<
215  "<KVElasticScatter::CalculateScattering> : Set projectile properties first"
216  << endl;
217  return;
218  }
219  if (!fEnergy) {
220  cout <<
221  "<KVElasticScatter::CalculateScattering> : Set projectile energy first"
222  << endl;
223  return;
224  }
225  if (!fDetector) {
226  cout <<
227  "<KVElasticScatter::CalculateScattering> : Set detector first" <<
228  endl;
229  return;
230  }
231  if (!fTarget) {
232  cout <<
233  "<KVElasticScatter::CalculateScattering> : No target set. Set run first."
234  << endl;
235  return;
236  }
237 
238  /* -------------------------------------------------------------------------------------------------------------------------- */
239 
240  //set up histograms
241 
242  /* -------------------------------------------------------------------------------------------------------------------------- */
243 
244  // -- histograms for debugging: depth in target and angular distribution
245  if (fDepth)
246  delete fDepth;
247  if (fTheta)
248  delete fTheta;
249  fDepth =
250  new TH1F("hDepth", "Depth (mg/cm2)", 500, 0.,
251  fTarget->GetTotalEffectiveThickness());
252  fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
253 
254  /* -------------------------------------------------------------------------------------------------------------------------- */
255 
256  //set up histograms for all detectors particle passes through
257  if (!fHistos) {
258  fHistos = new TObjArray(fAlignedDetectors.GetSize());
259  fHistos->SetOwner(); //will delete histos it stores
260  }
261  else {
262  fHistos->Clear(); //delete any previous histograms
263  }
264  KVDetector* d;
265  TIter n(&fAlignedDetectors);
266  while ((d = (KVDetector*) n())) {
267  fHistos->
268  Add(new
269  TH1F(Form("hEloss_%s", d->GetName()), "Eloss (MeV)", fBinE, 0., 0.));
270  }
271 
272  /* -------------------------------------------------------------------------------------------------------------------------- */
273 
274  //set up kinematics
275  if (!fKinematics)
276  fKinematics = new KV2Body;
277  fProj.SetEnergy(fEnergy);
278  fProj.SetTheta(0);
279 
280  /* -------------------------------------------------------------------------------------------------------------------------- */
281 
282  //set random interaction point for scattering
283  if (fIntLayer) {
284  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
285  }
286  else {
287  fTarget->GetInteractionPoint(&fProj);
288  }
289 
290  /* -------------------------------------------------------------------------------------------------------------------------- */
291 
292  //get target nucleus properties from scattering layer
293  TVector3 IP = fTarget->GetInteractionPoint();
294  KVMaterial* targ_mat = fTarget->GetLayer(IP);
295  KVNucleus t;
296  t.SetZ((Int_t) targ_mat->GetZ());
297  t.SetA((Int_t) targ_mat->GetMass());
298  fKinematics->SetTarget(t);
299 
300  /* -------------------------------------------------------------------------------------------------------------------------- */
301 
302  //set excited state of target nucleus - in other words, dissipated energy for
303  //reaction due to excitation of target
304  // check order !!!
305  fKinematics->SetEDiss(fExx);
306 
307  /* -------------------------------------------------------------------------------------------------------------------------- */
308 
309  Double_t xsec;
310  for (int i = 0; i < N; i++) {
311  //calculate slowing of incoming projectile
312  fTarget->SetIncoming();
313  fTarget->DetectParticle(&fProj);
314  fKinematics->SetProjectile(fProj);
315  if (!has_outgoing) fOutgoing = fProj;
316 
317  fKinematics->SetOutgoing(fOutgoing);
318 // else fKinematics->SetOutgoing(fProj);
319  fKinematics->CalculateKinematics();
320  //set random direction of outgoing projectile
321 
322  double th, ph;
323  th = ph = 0.;
324  fDetector->GetRandomAngles(th, ph);
325 
326  //set energy of scattered nucleus
327  //WARNING: for inverse kinematics reactions, their are two energies for
328  //each angle below the maximum scattering angle.
329  //We only use the highest energy corresponding to the most forward CM angle.
330  Double_t e1, e2;
331  fKinematics->GetELab(3, th, 3, e1, e2);
332  fOutgoing.SetEnergy(TMath::Max(e1, e2));
333  fOutgoing.SetTheta(th);
334  fOutgoing.SetPhi(ph);
335  xsec = TMath::Abs(fKinematics->GetXSecRuthLab(fOutgoing.GetTheta()));
336  fTheta->Fill(fOutgoing.GetTheta(), xsec);
337  //slowing of outgoing projectile in target
338  fTarget->SetOutgoing();
339 // fTarget->DetectParticle(&fProj);
340  fTarget->DetectParticle(&fOutgoing);
341  fDepth->Fill(IP.z());
342 
343  //now detect particle in detector(s)
344  int j = 0;
345  n.Reset();
346  while ((d = (KVDetector*) n())) {
347  auto cew = d->GetCentreOfEntranceWindow();
348  d->DetectParticle(&fOutgoing, &cew);
349  //fill histograms
350  ((TH1F*)(*fHistos)[j++])->Fill(d->GetEnergy(), xsec);
351  //prepare for next round: set energy loss to zero
352  d->Clear();
353  }
354  fProj.SetEnergy(fEnergy);
355  fProj.SetTheta(0);
356  fProj.GetParameters()->Clear();
357  fOutgoing.GetParameters()->Clear();
358  //set random interaction point for scattering
359  if (fIntLayer) {
360  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
361  }
362  else {
363  fTarget->GetInteractionPoint(&fProj);
364  //if target is multilayer and the interaction layer is not fixed,
365  //the layer & hence the target nucleus may change
366  if (fMultiLayer) {
367  targ_mat = fTarget->GetLayer(fTarget->GetInteractionPoint());
368  KVNucleus t;
369  t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
370  fKinematics->SetTarget(t);
371  }
372  }
373  IP = fTarget->GetInteractionPoint();
374  }
375 }
376 
377 
378 
379 
380 
383 
385 {
386  //Energy loss in detector with given 'label' through which scattered particle passes.
387 
388  return (fDetInd.HasParameter(label) ? GetEnergy(fDetInd.GetIntValue(label)) : 0);
389 }
390 
391 
392 
393 
400 
402 {
403  //Energy loss in any detector through which scattered particle passes.
404  //
405  //The index corresponds to the order in which detectors are crossed by the
406  //particle, beginning with 0 for the first detector, and ending with
407  //(GetNDets()-1)
408  return (TH1F*)(*fHistos)[index];
409 }
410 
411 
412 //__________________________________________________________________//
413 
int Int_t
#define d(i)
#define e(i)
char Char_t
constexpr Bool_t kFALSE
double Double_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]
char * Form(const char *fmt,...)
Relativistic binary kinematics calculator.
Definition: KV2Body.h:166
Calculate elastic scattering spectra in multidetector arrays.
void SetDetector(const Char_t *det)
TH1F * GetEnergy()
Return pointer to energy loss histogram for chosen detector (in MeV)
void CalculateScattering(Int_t N)
void SetEbinning(Int_t nbins=500)
void SetProjectile(Int_t Z, Int_t A)
Set projectile Z and A.
void SetRun(Int_t run)
Set detector parameters, target, etc. for run.
void SetTargetScatteringLayer(const Char_t *name)
void SetEnergy(Double_t E)
Set energy of projectile in MeV.
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:89
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
void SetA(Int_t a)
Definition: KVNucleus.cpp:647
void SetZ(Int_t z, Char_t mt=-1)
Definition: KVNucleus.cpp:696
void SetZandA(Int_t z, Int_t a)
Set atomic number and mass number.
Definition: KVNucleus.cpp:713
void SetRandomized(Bool_t r=kTRUE)
Definition: KVTarget.h:253
Double_t z() const
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,...)
Add
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
ClassImp(TPyArg)