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 
119 
121 {
122  //Set name of detector which will detect particle
123  fDetector = gMultiDetArray->GetDetector(det);
124  if (!fDetector) {
125  Error("SetDetector", "Detector %s is unknown!", det);
126  return;
127  }
128  //get list of all detectors particle must pass through to get to required detector
129  if (!fDetector->GetNode()->GetForwardTrajectories()) {
130  Error("SetDetector", "There is no forwards trajectory from detector %s:", det);
131  return;
132  }
133  if (fDetector->GetNode()->GetForwardTrajectories()->GetEntries() > 1) {
134  Error("SetDetector", "Ambiguous: there is more than one forwards trajectory from detector %s:", det);
135  fDetector->GetNode()->GetForwardTrajectories()->ls();
136  return;
137  }
138  auto traj = (KVGeoDNTrajectory*)fDetector->GetNode()->GetForwardTrajectories()->First();
139  //we store the association between detector type and index in list
140  fDetInd.Clear();
141  fAlignedDetectors.Clear();
142 
143  traj->IterateBackFrom();
144 
145  Int_t i = 0;
146  KVGeoDetectorNode* gdn = nullptr;
147  while ((gdn = traj->GetNextNode()) && gdn != fDetector->GetNode()) {
148  KVDetector* d = gdn->GetDetector();
149  fAlignedDetectors.Add(d);
150  //check same type is not already present
151  if (fDetInd.HasParameter(d->GetType())) {
152  //detetors with same type will be called "Type_1", "Type_2" etc.
153  TString newname;
154  int j = 1;
155  newname.Form("%s_%d", d->GetType(), j++);
156  while (fDetInd.HasParameter(newname.Data()))
157  newname.Form("%s_%d", d->GetType(), j++);
158  fDetInd.SetValue(newname.Data(), i);
159  }
160  else {
161  fDetInd.SetValue(d->GetType(), i);
162  }
163  ++i;
164  }
165  //store number of aligned detectors
166  fNDets = i;
167 }
168 
169 
170 
171 
177 
179 {
180  //For multilayer targets, use this method to choose in which
181  //layer the scattering will take place.
182  //If name="", reset any previous choice so that scattering
183  //can take place in any layer
184  if (!fTarget) {
185  cout <<
186  "<KVElasticScatter::SetTargetScatteringLayer> : No target set. Set run first."
187  << endl;
188  return;
189  }
190  fIntLayer = fTarget->GetLayerIndex(name);
191  if (fIntLayer)
192  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
193 }
194 
195 
196 
197 
201 
203 {
204  //Set binning of the GetEnergy histogram
205  // Default value is 500
206  fBinE = nbins;
207 }
208 
209 
210 
214 
216 {
217  //Perform scattering 'N' times for current values
218  //of particle Z, A and energy, target excited state, and detector.
219 
220  if (!fProj.IsDefined()) {
221  cout <<
222  "<KVElasticScatter::CalculateScattering> : Set projectile properties first"
223  << endl;
224  return;
225  }
226  if (!fEnergy) {
227  cout <<
228  "<KVElasticScatter::CalculateScattering> : Set projectile energy first"
229  << endl;
230  return;
231  }
232  if (!fDetector) {
233  cout <<
234  "<KVElasticScatter::CalculateScattering> : Set detector first" <<
235  endl;
236  return;
237  }
238  if (!fTarget) {
239  cout <<
240  "<KVElasticScatter::CalculateScattering> : No target set. Set run first."
241  << endl;
242  return;
243  }
244 
245  /* -------------------------------------------------------------------------------------------------------------------------- */
246 
247  //set up histograms
248 
249  /* -------------------------------------------------------------------------------------------------------------------------- */
250 
251  // -- histograms for debugging: depth in target and angular distribution
252  if (fDepth)
253  delete fDepth;
254  if (fTheta)
255  delete fTheta;
256  fDepth =
257  new TH1F("hDepth", "Depth (mg/cm2)", 500, 0.,
258  fTarget->GetTotalEffectiveThickness());
259  fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
260 
261  /* -------------------------------------------------------------------------------------------------------------------------- */
262 
263  //set up histograms for all detectors particle passes through
264  if (!fHistos) {
265  fHistos = new TObjArray(fAlignedDetectors.GetSize());
266  fHistos->SetOwner(); //will delete histos it stores
267  }
268  else {
269  fHistos->Clear(); //delete any previous histograms
270  }
271  KVDetector* d;
272  TIter n(&fAlignedDetectors);
273  while ((d = (KVDetector*) n())) {
274  fHistos->
275  Add(new
276  TH1F(Form("hEloss_%s", d->GetName()), "Eloss (MeV)", fBinE, 0., 0.));
277  }
278 
279  /* -------------------------------------------------------------------------------------------------------------------------- */
280 
281  //set up kinematics
282  if (!fKinematics)
283  fKinematics = new KV2Body;
284  fProj.SetEnergy(fEnergy);
285  fProj.SetTheta(0);
286 
287  /* -------------------------------------------------------------------------------------------------------------------------- */
288 
289  //set random interaction point for scattering
290  if (fIntLayer) {
291  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
292  }
293  else {
294  fTarget->GetInteractionPoint(&fProj);
295  }
296 
297  /* -------------------------------------------------------------------------------------------------------------------------- */
298 
299  //get target nucleus properties from scattering layer
300  TVector3 IP = fTarget->GetInteractionPoint();
301  KVMaterial* targ_mat = fTarget->GetLayer(IP);
302  KVNucleus t;
303  t.SetZ((Int_t) targ_mat->GetZ());
304  t.SetA((Int_t) targ_mat->GetMass());
305  fKinematics->SetTarget(t);
306 
307  /* -------------------------------------------------------------------------------------------------------------------------- */
308 
309  //set excited state of target nucleus - in other words, dissipated energy for
310  //reaction due to excitation of target
311  // check order !!!
312  fKinematics->SetEDiss(fExx);
313 
314  /* -------------------------------------------------------------------------------------------------------------------------- */
315 
316  Double_t xsec;
317  for (int i = 0; i < N; i++) {
318  //calculate slowing of incoming projectile
319  fTarget->SetIncoming();
320  fTarget->DetectParticle(&fProj);
321  fKinematics->SetProjectile(fProj);
322  if (!has_outgoing) fOutgoing = fProj;
323 
324  fKinematics->SetOutgoing(fOutgoing);
325 // else fKinematics->SetOutgoing(fProj);
326  fKinematics->CalculateKinematics();
327  //set random direction of outgoing projectile
328 
329  double th, ph;
330  th = ph = 0.;
331  fDetector->GetRandomAngles(th, ph);
332 
333  //set energy of scattered nucleus
334  //WARNING: for inverse kinematics reactions, their are two energies for
335  //each angle below the maximum scattering angle.
336  //We only use the highest energy corresponding to the most forward CM angle.
337  Double_t e1, e2;
338  fKinematics->GetELab(3, th, 3, e1, e2);
339  fOutgoing.SetEnergy(TMath::Max(e1, e2));
340  fOutgoing.SetTheta(th);
341  fOutgoing.SetPhi(ph);
342  xsec = TMath::Abs(fKinematics->GetXSecRuthLab(fOutgoing.GetTheta()));
343  fTheta->Fill(fOutgoing.GetTheta(), xsec);
344  //slowing of outgoing projectile in target
345  fTarget->SetOutgoing();
346 // fTarget->DetectParticle(&fProj);
347  fTarget->DetectParticle(&fOutgoing);
348  fDepth->Fill(IP.z());
349 
350  //now detect particle in detector(s)
351  int j = 0;
352  n.Reset();
353  while ((d = (KVDetector*) n())) {
354  auto cew = d->GetCentreOfEntranceWindow();
355  d->DetectParticle(&fOutgoing, &cew);
356  //fill histograms
357  ((TH1F*)(*fHistos)[j++])->Fill(d->GetEnergy(), xsec);
358  //prepare for next round: set energy loss to zero
359  d->Clear();
360  }
361  fProj.SetEnergy(fEnergy);
362  fProj.SetTheta(0);
363  fProj.GetParameters()->Clear();
364  fOutgoing.GetParameters()->Clear();
365  //set random interaction point for scattering
366  if (fIntLayer) {
367  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
368  }
369  else {
370  fTarget->GetInteractionPoint(&fProj);
371  //if target is multilayer and the interaction layer is not fixed,
372  //the layer & hence the target nucleus may change
373  if (fMultiLayer) {
374  targ_mat = fTarget->GetLayer(fTarget->GetInteractionPoint());
375  KVNucleus t;
376  t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
377  fKinematics->SetTarget(t);
378  }
379  }
380  IP = fTarget->GetInteractionPoint();
381  }
382 }
383 
384 
385 
386 
387 
393 
395 {
396  //Energy loss in detector of given 'type' through which scattered particle passes.
397  //Warning: if there are several detectors of the same type in the list of detectors
398  //through which the particle passes, the first one (as seen by the impinging
399  //particle) will have type "type", the second "type_1", the third "type_2", etc.
400 
401  return (fDetInd.HasParameter(type) ? GetEnergy(fDetInd.GetIntValue(type)) : 0);
402 }
403 
404 
405 
406 
412 
414 {
415  //Energy loss in any detector through which scattered particle passes.
416  //The index corresponds to the order in which detectors are crossed by the
417  //particle, beginning with 0 for the first detector, and ending with
418  //(GetNDets()-1)
419  return (TH1F*)(*fHistos)[index];
420 }
421 
422 
423 //__________________________________________________________________//
424 
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
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
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)
Set name of detector which will detect particle.
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
const char * Data() const
void Form(const char *fmt,...)
Double_t z() const
const Int_t n
void Error(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)