KaliVeda  1.13/01
Heavy-Ion Analysis Toolkit
KVElasticCountRates.cpp
Go to the documentation of this file.
1 /*
2 $Id: KVElasticCountRates.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 Fri Nov 20 2015
8 //Author: John Frankland
9 
10 #include "KVElasticCountRates.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 "KVNucleus.h"
19 #include "TObjArray.h"
20 #include "KVDatime.h"
21 
22 #include <TH2F.h>
23 #include <vector>
24 #include <algorithm>
25 
26 using namespace std;
27 
30 
31 
32 
35 KVElasticCountRates::KVElasticCountRates(Double_t theta_min, Double_t theta_max, Double_t phi_min, Double_t phi_max):
36  fAngularRange(theta_min, theta_max, phi_min, phi_max),
37  fBeamDirection(0, 0, 1)
38 {
39  //Default constructor
40  fDepth = fTheta = 0;
41  fBinE = 500;
42  fEnergy = 0;
43  fKinematics = 0;
44  fTarget = 0;
45  fIntLayer = 0;
46  fMultiLayer = kFALSE;
47  fVolume = (theta_max - theta_min) * (phi_max - phi_min) * TMath::DegToRad() * TMath::DegToRad();
48  //fVolume = (cos(theta_min*TMath::DegToRad())-cos(theta_max*TMath::DegToRad()))*(phi_max-phi_min)*TMath::DegToRad();
49  fExx = 0.;
50  if (gMultiDetArray) {
54  }
55  else {
56  Warning("KVElasticCountRates", "gMultiDetArray does not refer to a valid multidetector array");
57  printf("Define it before using this class, and put it in simulation mode : gMultiDetArray->SetSimMode(kTRUE)");
58  }
59 }
60 
61 
62 
63 
66 
67 KVElasticCountRates::~KVElasticCountRates()
68 {
69  //Destructor
70  if (fDepth)
71  delete fDepth;
72  if (fTheta)
73  delete fTheta;
74  if (fKinematics)
75  delete fKinematics;
76 }
77 
78 
79 
80 
83 
85 {
86  //Set detector parameters, target, etc. for run
90  fTarget = gMultiDetArray->GetTarget();
91  fAtomicDensity = fTarget->GetAtomsPerCM2() * 1.e-24; //in barn^-1 units
92  fTarget->SetRandomized();
93  fIntLayer = 0;
94  fMultiLayer = (fTarget->NumberOfLayers() > 1);
95 }
96 
97 
98 
101 
103 {
104  // Set projectile and beam energy [MeV/nucleon]
105  fProj = KVNucleus(nuc, e_sur_a);
106  fEnergy = fProj.GetE();
107 }
108 
109 
110 
111 
117 
119 {
120  //For multilayer targets, use this method to choose in which
121  //layer the scattering will take place.
122  //If name="", reset any previous choice so that scattering
123  //can take place in any layer
124  if (!fTarget) {
125  cout <<
126  "<KVElasticCountRates::SetTargetScatteringLayer> : No target set. Set run first."
127  << endl;
128  return;
129  }
130  fIntLayer = fTarget->GetLayerIndex(name);
131  if (fIntLayer)
132  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
133 }
134 
135 
136 
137 
142 
144 {
145  //Set the number of bins of the GetEnergy() histogram
146  //Default value is 500; this function has to be called before
147  //using CalculateScattering.
148  fBinE = nbins;
149 }
150 
151 
152 
158 
160 {
161  //Perform scattering 'N' times for current values
162  //of particle Z, A and energy, target excited state, and detector.
163  //Print out for each hit detector the mean energy loss & counting rate
164  //for a beam intensity of 10**7 pps
165 
166  fNtirages = N;
167 
168  fHistos.Delete();
169  fRates.clear();
170 
171  if (!fProj.IsDefined()) {
172  cout <<
173  "<KVElasticCountRates::CalculateScattering> : Set projectile properties first"
174  << endl;
175  return;
176  }
177  if (!fEnergy) {
178  cout <<
179  "<KVElasticCountRates::CalculateScattering> : Set projectile energy first"
180  << endl;
181  return;
182  }
183  if (!fTarget) {
184  cout <<
185  "<KVElasticCountRates::CalculateScattering> : No target set. Set run first."
186  << endl;
187  return;
188  }
189 
190  /* -------------------------------------------------------------------------------------------------------------------------- */
191 
192  //set up histograms
193 
194  /* -------------------------------------------------------------------------------------------------------------------------- */
195 
196  // -- histograms for debugging: depth in target and angular distribution
197  if (fDepth)
198  delete fDepth;
199  if (fTheta)
200  delete fTheta;
201  fDepth =
202  new TH1F("hDepth", "Depth (mg/cm2)", 500, 0., fTarget->GetTotalEffectiveThickness());
203  fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
204 
205  /* -------------------------------------------------------------------------------------------------------------------------- */
206 
207  //set up kinematics
208  if (!fKinematics)
209  fKinematics = new KV2Body;
210  fProj.SetEnergy(fEnergy);
211  fProj.SetTheta(0);
212 
213  /* -------------------------------------------------------------------------------------------------------------------------- */
214 
215  //set random interaction point for scattering
216  if (fIntLayer) {
217  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
218  }
219  else {
220  fTarget->GetInteractionPoint(&fProj);
221  }
222 
223  /* -------------------------------------------------------------------------------------------------------------------------- */
224 
225  //get target nucleus properties from scattering layer
226  TVector3 IP = fTarget->GetInteractionPoint();
227  KVMaterial* targ_mat = fTarget->GetLayer(IP);
228  KVNucleus t;
229  t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
230  fKinematics->SetTarget(t);
231 
232  /* -------------------------------------------------------------------------------------------------------------------------- */
233 
234  //set excited state of target nucleus - in other words, dissipated energy for
235  //reaction due to excitation of target
236  fKinematics->SetEDiss(fExx);
237 
238  /* -------------------------------------------------------------------------------------------------------------------------- */
239 
240  for (int i = 0; i < N; i++) {
241  //calculate slowing of incoming projectile
242  fTarget->SetIncoming();
243  fTarget->DetectParticle(&fProj);
244  fKinematics->SetProjectile(fProj);
245  fKinematics->SetOutgoing(fProj);
246  fKinematics->CalculateKinematics();
247  //set random direction of outgoing projectile
248  fAngularRange.GetRandomAngles(theta, phi, "random");
249  xsec = TMath::Abs(fKinematics->GetXSecRuthLab(theta));
250  //set energy of scattered nucleus
251  //WARNING: for inverse kinematics reactions, their are two energies for
252  //each angle below the maximum scattering angle.
253  //We only use the highest energy corresponding to the most forward CM angle.
254  Double_t e1, e2;
255  fKinematics->GetELab(3, theta, 3, e1, e2);
256  fProj.SetEnergy(TMath::Max(e1, e2));
257  fProj.SetTheta(theta);
258  fProj.SetPhi(phi);
259  fTheta->Fill(theta, xsec);
260  //slowing of outgoing projectile in target
261  fTarget->SetOutgoing();
262  fTarget->DetectParticle(&fProj);
263  //now detect particle in array
264  KVNameValueList* detectors = gMultiDetArray->DetectParticle(&fProj);
265  //fill histograms
266  fDepth->Fill(IP.z());
267  FillHistograms(detectors);
268  fProj.GetParameters()->Clear();
269  fProj.SetEnergy(fEnergy);
270  fProj.SetTheta(0);
271  //set random interaction point for scattering
272  if (fIntLayer) {
273  fTarget->SetInteractionLayer(fIntLayer, fBeamDirection);
274  }
275  else {
276  fTarget->GetInteractionPoint(&fProj);
277  //if target is multilayer and the interaction layer is not fixed,
278  //the layer & hence the target nucleus may change
279  if (fMultiLayer) {
280  targ_mat = fTarget->GetLayer(fTarget->GetInteractionPoint());
281  KVNucleus t;
282  t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
283  fKinematics->SetTarget(t);
284  }
285  }
286  IP = fTarget->GetInteractionPoint();
287  }
288 
289  PrintResults();
290 }
291 
292 
293 
294 
295 
296 
302 
304 {
305  // parse the list dets
306  // fill histograms with energy loss for all detectors
307  // clear the detector energy losses
308  // delete the list
309 
310  if (!dets) return;
311 
312  Int_t ndets = dets->GetNpar();
313  for (int i = 0; i < ndets; i++) {
314  TString detname = dets->GetNameAt(i);
315  KVDetector* det = gMultiDetArray->GetDetector(detname);
316  if (!det) continue;
317  TH1F* histo = (TH1F*)fHistos.FindObject(detname);
318  if (!histo) {
319  histo = new TH1F(detname, Form("Eloss in %s", detname.Data()), fBinE, 0, 0);
320  fHistos.Add(histo);
321  }
322  double de = dets->GetDoubleValue(i);
323  histo->Fill(de, xsec * sin(theta * TMath::DegToRad()));
324  histo = (TH1F*)fHistos.FindObject(detname + "_dW");
325  if (!histo) {
326  histo = new TH1F(detname + "_dW", Form("Solid angle of %s", detname.Data()), fBinE, 0, 0);
327  fHistos.Add(histo);
328  }
329  histo->Fill(de, sin(theta * TMath::DegToRad()));
330  TH2F* histo2 = (TH2F*)fHistos.FindObject(detname + "_map");
331  if (!histo2) {
332  histo2 = new TH2F(detname + "_map", Form("Map of %s", detname.Data()), 100, 0, 0, 100, 0, 0);
333  fHistos.Add(histo2);
334  }
335  histo2->Fill(theta, phi, xsec);
336  det->Clear();
337  }
338  delete dets;
339 }
340 
341 
345 
350 struct count_rate {
352  double counts;
353  double energy;
354  double theta;
355  double phi;
356  double fluence;
357  double dissipation;
358  double intXsec;
359  count_rate(TString n, double c, double e, double t, double p, double f, double d, double i)
360  : detector(n), counts(c), energy(e), theta(t), phi(p), fluence(f), dissipation(d), intXsec(i) {}
361  void print()
362  {
363  printf("%s \t: N=%8.2f/sec. \t <E>=%7.1f MeV \t Tot.Xsec=%9.3E barn \t fluence=%9.3E/sec./cm**2 \t dissip.=%9.3E MeV/sec./cm**2\n",
364  detector.Data(), counts, energy, intXsec, fluence, dissipation);
365  }
367  {
368  if (nl.HasParameter("DetList")) {
369  KVString tmp = nl.GetStringValue("DetList");
370  tmp += ",";
371  tmp += detector;
372  nl.SetValue("DetList", tmp.Data());
373  }
374  else {
375  nl.SetValue("DetList", detector.Data());
376  }
377  nl.SetValue(Form("%s.Counts(/sec)", detector.Data()), Form("%8.2f", counts));
378 
379  nl.SetValue(Form("%s.Counts(/sec)", detector.Data()), Form("%8.2f", counts));
380  nl.SetValue(Form("%s.Emean(MeV)", detector.Data()), Form("%7.1f", energy));
381  nl.SetValue(Form("%s.TotXsec(barn)", detector.Data()), Form("%9.3E", intXsec));
382  nl.SetValue(Form("%s.Fluence(/sec./cm**2)", detector.Data()), Form("%9.3E", fluence));
383  nl.SetValue(Form("%s.dissip(MeV/sec./cm**2)", detector.Data()), Form("%9.3E", dissipation));
384  }
385 };
386 
387 
389 
391 {
392  if (TMath::Abs(a.theta - b.theta) < 0.5) return a.phi < b.phi;
393  return a.theta < b.theta;
394 }
395 
396 
397 
400 
402 {
403  // Print mean energy deposit & counting rate for given beam intensity in particles per second
404 
405  TIter it(&fHistos);
406  TH1F* h;
407  fRates.clear();
408 
409  std::vector<count_rate> count_rates;
410 
411  while ((h = (TH1F*)it())) {
412  TString name = h->GetName();
413  if (!name.EndsWith("_dW") && !name.EndsWith("_map")) {
414  TH2F* map = (TH2F*)fHistos.FindObject(name + "_map");
415  // integrated cross-section
416  double intXsec = h->Integral() * fVolume / fNtirages;
417  // counting rate
418  double rate = fAtomicDensity * beam_intensity * intXsec;
419  // mean energy
420  double emean = h->GetMean();
421  KVDetector* det = gMultiDetArray->GetDetector(name);
422  double fluence = rate / det->GetEntranceWindowSurfaceArea();
423  double dissipation = emean * rate / det->GetEntranceWindowSurfaceArea();
424  count_rates.push_back(
425  count_rate(name, rate, emean, map->GetMean(), map->GetMean(2), fluence, dissipation, intXsec)
426  );
427  fRates[name.Data()] = KVElasticCountRate(rate, emean, intXsec, fluence, dissipation);
428  }
429  }
430  std::sort(count_rates.begin(), count_rates.end(), compare_count_rates);
431 
432  for (std::vector<count_rate>::iterator it = count_rates.begin(); it != count_rates.end(); ++it) {
433  it->print();
434  }
435 }
436 
437 
438 
441 
443 {
444  // Print mean energy deposit & counting rate for given beam intensity in particles per second
445  KVNameValueList nl;
446  nl.SetName("Generated by KVElasticCountRates::PutResultsInList method");
447  KVDatime dt;
448  nl.SetTitle(dt.AsSQLString());
449  TIter it(&fHistos);
450  TH1F* h;
451  fRates.clear();
452  nl.SetValue("Intensity(pps)", beam_intensity);
453 
454  std::vector<count_rate> count_rates;
455 
456  while ((h = (TH1F*)it())) {
457  TString name = h->GetName();
458  if (!name.EndsWith("_dW") && !name.EndsWith("_map")) {
459  TH2F* map = (TH2F*)fHistos.FindObject(name + "_map");
460  // integrated cross-section
461  double intXsec = h->Integral() * fVolume / fNtirages;
462  // counting rate
463  double rate = fAtomicDensity * beam_intensity * intXsec;
464  // mean energy
465  double emean = h->GetMean();
466  KVDetector* det = gMultiDetArray->GetDetector(name);
467  double fluence = rate / det->GetEntranceWindowSurfaceArea();
468  double dissipation = emean * rate / det->GetEntranceWindowSurfaceArea();
469  count_rates.push_back(
470  count_rate(name, rate, emean, map->GetMean(), map->GetMean(2), fluence, dissipation, intXsec)
471  );
472  fRates[name.Data()] = KVElasticCountRate(rate, emean, intXsec, fluence, dissipation);
473  }
474  }
475  std::sort(count_rates.begin(), count_rates.end(), compare_count_rates);
476 
477  for (std::vector<count_rate>::iterator it = count_rates.begin(); it != count_rates.end(); ++it) {
478  it->PutInList(nl);
479  }
480  return nl;
481 }
482 
483 
484 //__________________________________________________________________//
485 
int Int_t
bool compare_count_rates(count_rate a, count_rate b)
KVMultiDetArray * gMultiDetArray
ClassImp(KVPartitionList) void KVPartitionList
Initialisation.
static Double_t energy[]
#define d(i)
#define b(i)
#define f(i)
#define c(i)
#define e(i)
char Char_t
const Bool_t kFALSE
double Double_t
const Bool_t kTRUE
#define N
double sin(double)
char * Form(const char *fmt,...)
Relativistic binary kinematics calculator.
Definition: KV2Body.h:165
Extension of TDatime to handle various useful date formats.
Definition: KVDatime.h:32
Base class for detector geometry description.
Definition: KVDetector.h:159
virtual void Clear(Option_t *opt="")
Definition: KVDetector.cpp:596
virtual Double_t GetEntranceWindowSurfaceArea()
Return surface area of first layer of detector in cm2.
Calculate elastic scattering count rates in multidetector arrays.
void PrintResults(Double_t beam_intensity=1.e+07)
Print mean energy deposit & counting rate for given beam intensity in particles per second.
KVNameValueList PutResultsInList(Double_t beam_intensity=1.e+07)
Print mean energy deposit & counting rate for given beam intensity in particles per second.
void FillHistograms(KVNameValueList *)
void SetProjectile(const Char_t *nuc, Double_t e_sur_a)
Set projectile and beam energy [MeV/nucleon].
void SetEbinning(Int_t nbins=500)
void CalculateScattering(Int_t N=10000)
void SetRun(Int_t run)
Set detector parameters, target, etc. for run.
void SetTargetScatteringLayer(const Char_t *name)
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:93
Double_t GetZ() const
Definition: KVMaterial.cpp:393
Double_t GetMass() const
Definition: KVMaterial.cpp:305
virtual KVNameValueList * DetectParticle(KVNucleus *part)
KVTarget * GetTarget()
virtual void SetParameters(UInt_t n, Bool_t physics_parameters_only=kFALSE)
void SetFilterType(Int_t t)
virtual void InitializeIDTelescopes()
virtual void SetROOTGeometry(Bool_t on=kTRUE)
virtual void SetSimMode(Bool_t on=kTRUE)
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Double_t GetDoubleValue(const Char_t *name) const
void SetValue(const Char_t *name, value_type value)
const Char_t * GetNameAt(Int_t idx) const
Int_t GetNpar() const
return the number of stored parameters
const Char_t * GetStringValue(const Char_t *name) const
Bool_t HasParameter(const Char_t *name) const
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:125
void SetZandA(Int_t z, Int_t a)
Set atomic number and mass number.
Definition: KVNucleus.cpp:721
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition: KVString.h:72
Double_t GetAtomsPerCM2() const
virtual UInt_t GetUnits() const;
Definition: KVTarget.cpp:926
const char * AsSQLString() const
virtual Double_t GetMean(Int_t axis=1) const
virtual TObject * FindObject(const char *name) const
virtual Int_t Fill(const char *name, Double_t w)
virtual Double_t Integral(Int_t binx1, Int_t binx2, Option_t *option="") const
virtual Int_t Fill(const char *namex, const char *namey, Double_t w)
virtual const char * GetName() const
virtual void SetTitle(const char *title="")
virtual void SetName(const char *name)
Bool_t EndsWith(const char *pat, ECaseCompare cmp=kExact) const
const char * Data() const
Double_t z() const
const Int_t n
TH1 * h
void Warning(const char *location, const char *va_(fmt),...)
constexpr Double_t DegToRad()
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
Utility class used by KVElasticCountRates to store results.
Utility class used by KVElasticCountRates.
count_rate(TString n, double c, double e, double t, double p, double f, double d, double i)
void PutInList(KVNameValueList &nl)
auto * a