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  fAlignedDetectors = 0;
40  fExx = 0.;
41  fHistos = 0;
42  fDetInd = 0;
43  if (gMultiDetArray) {
44  gMultiDetArray->SetSimMode(kTRUE);
45  }
46  else {
47  Warning("KVElasticScatter", "gMultiDetArray does not refer to a valid multidetector array");
48  printf("Define it before using this class, and put it in simulation mode : gMultiDetArray->SetSimMode(kTRUE)");
49  }
50 }
51 
52 
53 
54 
57 
58 KVElasticScatter::~KVElasticScatter()
59 {
60  //Destructor
61  if (fDepth)
62  delete fDepth;
63  if (fTheta)
64  delete fTheta;
65  if (fKinematics)
66  delete fKinematics;
67  if (fHistos)
68  delete fHistos;
69  if (fDetInd)
70  delete fDetInd;
71  gMultiDetArray->SetSimMode(kFALSE);
72 }
73 
74 
75 
76 
79 
81 {
82  //Set detector parameters, target, etc. for run
83  gMultiDetArray->SetParameters(run);
84  fTarget = gMultiDetArray->GetTarget();
85  fTarget->SetRandomized();
86  fIntLayer = 0;
87  fMultiLayer = (fTarget->NumberOfLayers() > 1);
88 }
89 
90 
91 
92 
95 
97 {
98  //Set projectile Z and A
99 
100  fProj.SetZandA(Z, A);
101 }
102 
103 
104 
105 
108 
110 {
111  //Set energy of projectile in MeV
112 
113  fProj.SetEnergy(e);
114  fEnergy = e;
115 }
116 
117 
118 
119 
122 
124 {
125  //Set name of detector which will detect particle
126  fDetector = gMultiDetArray->GetDetector(det);
127  //get list of all detectors particle must pass through to get to required detector
128  fAlignedDetectors = nullptr;
129  Fatal("SetDetector", "Reimplement using particle trajectories");
130 // fDetector->GetAlignedDetectors(KVGroup::kForwards);
131  //we store the association between detector type and index in list
132  if (!fDetInd)
133  fDetInd = new KVNameValueList;
134  else
135  fDetInd->Clear();
136  KVDetector* d;
137  TIter n(fAlignedDetectors);
138  Int_t i = 0;
139  while ((d = (KVDetector*) n())) {
140  //check same type is not already present
141  if (fDetInd->HasParameter(d->GetType())) {
142  //detetors with same type will be called "Type_1", "Type_2" etc.
143  TString newname;
144  int j = 1;
145  newname.Form("%s_%d", d->GetType(), j++);
146  while (fDetInd->HasParameter(newname.Data()))
147  newname.Form("%s_%d", d->GetType(), j++);
148  fDetInd->SetValue(newname.Data(), i);
149  }
150  else {
151  fDetInd->SetValue(d->GetType(), i);
152  }
153  i++;
154  }
155  //store number of aligned detectors
156  fNDets = i;
157 }
158 
159 
160 
161 
167 
169 {
170  //For multilayer targets, use this method to choose in which
171  //layer the scattering will take place.
172  //If name="", reset any previous choice so that scattering
173  //can take place in any layer
174  if (!fTarget) {
175  cout <<
176  "<KVElasticScatter::SetTargetScatteringLayer> : No target set. Set run first."
177  << endl;
178  return;
179  }
180  fIntLayer = fTarget->GetLayerIndex(name);
181  if (fIntLayer)
182  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
183 }
184 
185 
186 
187 
191 
193 {
194  //Set binning of the GetEnergy histogram
195  // Default value is 500
196  fBinE = nbins;
197 }
198 
199 
200 
204 
206 {
207  //Perform scattering 'N' times for current values
208  //of particle Z, A and energy, target excited state, and detector.
209 
210  if (!fProj.IsDefined()) {
211  cout <<
212  "<KVElasticScatter::CalculateScattering> : Set projectile properties first"
213  << endl;
214  return;
215  }
216  if (!fEnergy) {
217  cout <<
218  "<KVElasticScatter::CalculateScattering> : Set projectile energy first"
219  << endl;
220  return;
221  }
222  if (!fDetector) {
223  cout <<
224  "<KVElasticScatter::CalculateScattering> : Set detector first" <<
225  endl;
226  return;
227  }
228  if (!fTarget) {
229  cout <<
230  "<KVElasticScatter::CalculateScattering> : No target set. Set run first."
231  << endl;
232  return;
233  }
234 
235  /* -------------------------------------------------------------------------------------------------------------------------- */
236 
237  //set up histograms
238 
239  /* -------------------------------------------------------------------------------------------------------------------------- */
240 
241  // -- histograms for debugging: depth in target and angular distribution
242  if (fDepth)
243  delete fDepth;
244  if (fTheta)
245  delete fTheta;
246  fDepth =
247  new TH1F("hDepth", "Depth (mg/cm2)", 500, 0.,
248  fTarget->GetTotalEffectiveThickness());
249  fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
250 
251  /* -------------------------------------------------------------------------------------------------------------------------- */
252 
253  //set up histograms for all detectors particle passes through
254  if (!fHistos) {
255  fHistos = new TObjArray(fAlignedDetectors->GetSize());
256  fHistos->SetOwner(); //will delete histos it stores
257  }
258  else {
259  fHistos->Clear(); //delete any previous histograms
260  }
261  KVDetector* d;
262  TIter n(fAlignedDetectors);
263  while ((d = (KVDetector*) n())) {
264  fHistos->
265  Add(new
266  TH1F(Form("hEloss_%s", d->GetName()), "Eloss (MeV)", fBinE, 0., 0.));
267  }
268 
269  /* -------------------------------------------------------------------------------------------------------------------------- */
270 
271  //set up kinematics
272  if (!fKinematics)
273  fKinematics = new KV2Body;
274  fProj.SetEnergy(fEnergy);
275  fProj.SetTheta(0);
276 
277  /* -------------------------------------------------------------------------------------------------------------------------- */
278 
279  //set random interaction point for scattering
280  if (fIntLayer) {
281  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
282  }
283  else {
284  fTarget->GetInteractionPoint(&fProj);
285  }
286 
287  /* -------------------------------------------------------------------------------------------------------------------------- */
288 
289  //get target nucleus properties from scattering layer
290  TVector3 IP = fTarget->GetInteractionPoint();
291  KVMaterial* targ_mat = fTarget->GetLayer(IP);
292  KVNucleus t;
293  t.SetZ((Int_t) targ_mat->GetZ());
294  t.SetA((Int_t) targ_mat->GetMass());
295  fKinematics->SetTarget(t);
296 
297  /* -------------------------------------------------------------------------------------------------------------------------- */
298 
299  //set excited state of target nucleus - in other words, dissipated energy for
300  //reaction due to excitation of target
301  fKinematics->SetEDiss(fExx);
302 
303  /* -------------------------------------------------------------------------------------------------------------------------- */
304 
305  Double_t xsec;
306  for (int i = 0; i < N; i++) {
307  //calculate slowing of incoming projectile
308  fTarget->SetIncoming();
309  fTarget->DetectParticle(&fProj);
310  fKinematics->SetProjectile(fProj);
311  fKinematics->SetOutgoing(fProj);
312  fKinematics->CalculateKinematics();
313  //set random direction of outgoing projectile
314 
315  double th, ph;
316  th = ph = 0.;
317  fDetector->GetRandomAngles(th, ph);
318 
319  //set energy of scattered nucleus
320  //WARNING: for inverse kinematics reactions, their are two energies for
321  //each angle below the maximum scattering angle.
322  //We only use the highest energy corresponding to the most forward CM angle.
323  Double_t e1, e2;
324  fKinematics->GetELab(3, th, 3, e1, e2);
325  fProj.SetEnergy(TMath::Max(e1, e2));
326  fProj.SetTheta(th);
327  fProj.SetPhi(ph);
328  xsec = TMath::Abs(fKinematics->GetXSecRuthLab(fProj.GetTheta()));
329  fTheta->Fill(fProj.GetTheta(), xsec);
330  //slowing of outgoing projectile in target
331  fTarget->SetOutgoing();
332  fTarget->DetectParticle(&fProj);
333  //now detect particle in detector(s)
334  fAlignedDetectors->R__FOR_EACH(KVDetector, DetectParticle)(&fProj);
335  //fill histograms
336  fDepth->Fill(IP.z());
337  int j = 0;
338  n.Reset();
339  while ((d = (KVDetector*) n())) {
340  ((TH1F*)(*fHistos)[j++])->Fill(d->GetEnergy(), xsec);
341  //prepare for next round: set energy loss to zero
342  d->Clear();
343  }
344  fProj.SetEnergy(fEnergy);
345  fProj.SetTheta(0);
346  fProj.GetParameters()->Clear();
347  //set random interaction point for scattering
348  if (fIntLayer) {
349  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
350  }
351  else {
352  fTarget->GetInteractionPoint(&fProj);
353  //if target is multilayer and the interaction layer is not fixed,
354  //the layer & hence the target nucleus may change
355  if (fMultiLayer) {
356  targ_mat = fTarget->GetLayer(fTarget->GetInteractionPoint());
357  KVNucleus t;
358  t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
359  fKinematics->SetTarget(t);
360  }
361  }
362  IP = fTarget->GetInteractionPoint();
363  }
364 }
365 
366 
367 
368 
369 
375 
377 {
378  //Energy loss in detector of given 'type' through which scattered particle passes.
379  //Warning: if there are several detectors of the same type in the list of detectors
380  //through which the particle passes, the first one (as seen by the impinging
381  //particle) will have type "type", the second "type_1", the third "type_2", etc.
382 
383  return (fDetInd->HasParameter(type) ? GetEnergy(fDetInd->GetIntValue(type)) : 0);
384 }
385 
386 
387 
388 
394 
396 {
397  //Energy loss in any detector through which scattered particle passes.
398  //The index corresponds to the order in which detectors are crossed by the
399  //particle, beginning with 0 for the first detector, and ending with
400  //(GetNDets()-1)
401  return (TH1F*)(*fHistos)[index];
402 }
403 
404 
405 //__________________________________________________________________//
406 
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.
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)
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
void Clear(Option_t *opt="") override
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 Fatal(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)