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