KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVElasticCountRates.cpp
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
26using namespace std;
27
30
31
32
34
36 const TVector3& beam_dir):
37 fAngularRange(theta_min, theta_max, phi_min, phi_max),
38 fBeamDirection{beam_dir}
39{
40 //Default constructor
41 fDepth = fTheta = 0;
42 fBinE = 500;
43 fEnergy = 0;
44 fKinematics = 0;
45 fTarget = 0;
46 fIntLayer = 0;
47 fMultiLayer = kFALSE;
48 fVolume = (theta_max - theta_min) * (phi_max - phi_min) * TMath::DegToRad() * TMath::DegToRad();
49 //fVolume = (cos(theta_min*TMath::DegToRad())-cos(theta_max*TMath::DegToRad()))*(phi_max-phi_min)*TMath::DegToRad();
50 fExx = 0.;
51 if (gMultiDetArray) {
52 gMultiDetArray->SetSimMode(kTRUE);
54 gMultiDetArray->InitializeIDTelescopes();
55 }
56 else {
57 Warning("KVElasticCountRates", "gMultiDetArray does not refer to a valid multidetector array");
58 printf("Define it before using this class, and put it in simulation mode : gMultiDetArray->SetSimMode(kTRUE)");
59 }
60}
61
62
63
64
67
69{
70 //Destructor
71 if (fDepth)
72 delete fDepth;
73 if (fTheta)
74 delete fTheta;
75 if (fKinematics)
76 delete fKinematics;
77}
78
79
80
81
84
86{
87 //Set detector parameters, target, etc. for run
88 gMultiDetArray->SetParameters(run);
89 gMultiDetArray->InitializeIDTelescopes();
90 gMultiDetArray->SetROOTGeometry();
91 fTarget = gMultiDetArray->GetTarget();
92 fAtomicDensity = fTarget->GetAtomsPerCM2() * 1.e-24; //in barn^-1 units
94 fIntLayer = 0;
96}
97
98
99
102
104{
105 // Set projectile and beam energy [MeV/nucleon] & direction
106 fProj = KVNucleus(nuc, e_sur_a);
107 fEnergy = fProj.GetE();
109}
110
111
112
113
119
121{
122 //For multilayer targets, use this method to choose in which
123 //layer the scattering will take place.
124 //If name="", reset any previous choice so that scattering
125 //can take place in any layer
126 if (!fTarget) {
127 cout <<
128 "<KVElasticCountRates::SetTargetScatteringLayer> : No target set. Set run first."
129 << endl;
130 return;
131 }
133 if (fIntLayer)
135}
136
137
138
139
144
146{
147 //Set the number of bins of the GetEnergy() histogram
148 //Default value is 500; this function has to be called before
149 //using CalculateScattering.
150 fBinE = nbins;
151}
152
153
154
160
162{
163 //Perform scattering 'N' times for current values
164 //of particle Z, A and energy, target excited state, and detector.
165 //Print out for each hit detector the mean energy loss & counting rate
166 //for a beam intensity of 10**7 pps
167
168 fNtirages = N;
169
170 fHistos.Delete();
171 fRates.clear();
172
173 if (!fProj.IsDefined()) {
174 cout <<
175 "<KVElasticCountRates::CalculateScattering> : Set projectile properties first"
176 << endl;
177 return;
178 }
179 if (!fEnergy) {
180 cout <<
181 "<KVElasticCountRates::CalculateScattering> : Set projectile energy first"
182 << endl;
183 return;
184 }
185 if (!fTarget) {
186 cout <<
187 "<KVElasticCountRates::CalculateScattering> : No target set. Set run first."
188 << endl;
189 return;
190 }
191
192 /* -------------------------------------------------------------------------------------------------------------------------- */
193
194 //set up histograms
195
196 /* -------------------------------------------------------------------------------------------------------------------------- */
197
198 // -- histograms for debugging: depth in target and angular distribution
199 if (fDepth)
200 delete fDepth;
201 if (fTheta)
202 delete fTheta;
203 fDepth =
204 new TH1F("hDepth", "Depth (mg/cm2)", 500, 0., fTarget->GetTotalEffectiveThickness());
205 fTheta = new TH1F("hTheta", "Theta (deg.)", 500, 0., 0.);
206 fGlobalMap = new TH2F("hGlobalMap", "X-section in x(vert)-y(hori) plane",
209
210 /* -------------------------------------------------------------------------------------------------------------------------- */
211
212 //set up kinematics
213 if (!fKinematics)
214 fKinematics = new KV2Body;
215
216 /* -------------------------------------------------------------------------------------------------------------------------- */
217
218 //set random interaction point for scattering
219 if (fIntLayer) {
221 }
222 else {
224 }
225
226 /* -------------------------------------------------------------------------------------------------------------------------- */
227
228 //get target nucleus properties from scattering layer
230 KVMaterial* targ_mat = fTarget->GetLayer(IP);
231 KVNucleus t;
232 t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
234
235 /* -------------------------------------------------------------------------------------------------------------------------- */
236
237 //set excited state of target nucleus - in other words, dissipated energy for
238 //reaction due to excitation of target
240
241 /* -------------------------------------------------------------------------------------------------------------------------- */
242
243 for (int i = 0; i < N; i++) {
244 //calculate slowing of incoming projectile
250 //set random direction of outgoing projectile
254 // get Rutherford Xsec for projectile scattered to theta degrees in lab
255 // taking into account the initial direction of the beam (if beam direction = (0,0,1)
256 // this is just the same as theta)
257 double angle_with_beam = TMath::RadToDeg() * fProj.Angle(fBeamDirection);
258 xsec = TMath::Abs(fKinematics->GetXSecRuthLab(angle_with_beam));
259 //set energy of scattered nucleus
260 //WARNING: for inverse kinematics reactions, their are two energies for
261 //each angle below the maximum scattering angle.
262 //We only use the highest energy corresponding to the most forward CM angle.
263 Double_t e1, e2;
264 fKinematics->GetELab(3, angle_with_beam, 3, e1, e2);
265 fProj.SetEnergy(TMath::Max(e1, e2));
267 //slowing of outgoing projectile in target
270 //now detect particle in array
271 KVNameValueList* detectors = gMultiDetArray->DetectParticle(&fProj);
272 //fill histograms
273 fDepth->Fill(IP.z());
274 FillHistograms(detectors);
278 //set random interaction point for scattering
279 if (fIntLayer) {
281 }
282 else {
284 //if target is multilayer and the interaction layer is not fixed,
285 //the layer & hence the target nucleus may change
286 if (fMultiLayer) {
288 KVNucleus t;
289 t.SetZandA((Int_t) targ_mat->GetZ(), (Int_t) targ_mat->GetMass());
291 }
292 }
294 }
295
296 PrintResults();
297}
298
299
300
301
302
303
309
311{
312 // parse the list dets
313 // fill histograms with energy loss for all detectors
314 // clear the detector energy losses
315 // delete the list
316
317 if (!dets) return;
318
319 Int_t ndets = dets->GetNpar();
320 for (int i = 0; i < ndets; i++) {
321 TString detname = dets->GetNameAt(i);
322 KVDetector* det = gMultiDetArray->GetDetector(detname);
323 if (!det) continue;
324 TH1F* histo = (TH1F*)fHistos.FindObject(detname);
325 if (!histo) {
326 histo = new TH1F(detname, Form("Eloss in %s", detname.Data()), fBinE, 0, 0);
327 fHistos.Add(histo);
328 }
329 double de = dets->GetDoubleValue(i);
330 histo->Fill(de, xsec * sin(theta * TMath::DegToRad()));
331 histo = (TH1F*)fHistos.FindObject(detname + "_dW");
332 if (!histo) {
333 histo = new TH1F(detname + "_dW", Form("Solid angle of %s", detname.Data()), fBinE, 0, 0);
334 fHistos.Add(histo);
335 }
336 histo->Fill(de, sin(theta * TMath::DegToRad()));
337 TH2F* histo2 = (TH2F*)fHistos.FindObject(detname + "_map");
338 if (!histo2) {
339 histo2 = new TH2F(detname + "_map", Form("Map of %s", detname.Data()), 100, 0, 0, 100, 0, 0);
340 fHistos.Add(histo2);
341 }
342 histo2->Fill(theta, phi, xsec);
344 det->Clear();
345 }
346 delete dets;
347}
348
349
353
360 double counts;
361 double energy;
362 double theta;
363 double phi;
364 double fluence;
366 double intXsec;
367 count_rate(TString n, double c, double e, double t, double p, double f, double d, double i)
368 : detector(n), counts(c), energy(e), theta(t), phi(p), fluence(f), dissipation(d), intXsec(i) {}
369 void print()
370 {
371 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",
373 }
375 {
376 if (nl.HasParameter("DetList")) {
377 KVString tmp = nl.GetStringValue("DetList");
378 tmp += ",";
379 tmp += detector;
380 nl.SetValue("DetList", tmp.Data());
381 }
382 else {
383 nl.SetValue("DetList", detector.Data());
384 }
385 nl.SetValue(Form("%s.Counts(/sec)", detector.Data()), Form("%8.2f", counts));
386
387 nl.SetValue(Form("%s.Counts(/sec)", detector.Data()), Form("%8.2f", counts));
388 nl.SetValue(Form("%s.Emean(MeV)", detector.Data()), Form("%7.1f", energy));
389 nl.SetValue(Form("%s.TotXsec(barn)", detector.Data()), Form("%9.3E", intXsec));
390 nl.SetValue(Form("%s.Fluence(/sec./cm**2)", detector.Data()), Form("%9.3E", fluence));
391 nl.SetValue(Form("%s.dissip(MeV/sec./cm**2)", detector.Data()), Form("%9.3E", dissipation));
392 }
393};
394
395
397
398bool compare_count_rates(count_rate a, count_rate b)
399{
400 if (TMath::Abs(a.theta - b.theta) < 0.5) return a.phi < b.phi;
401 return a.theta < b.theta;
402}
403
404
405
408
410{
411 // Print mean energy deposit & counting rate for given beam intensity in particles per second
412
413 TIter it(&fHistos);
414 TH1F* h;
415 fRates.clear();
416
417 std::vector<count_rate> count_rates;
418
419 while ((h = (TH1F*)it())) {
420 TString name = h->GetName();
421 if (!name.EndsWith("_dW") && !name.EndsWith("_map")) {
422 TH2F* map = (TH2F*)fHistos.FindObject(name + "_map");
423 // integrated cross-section
424 double intXsec = h->Integral() * fVolume / fNtirages;
425 // counting rate
426 double rate = fAtomicDensity * beam_intensity * intXsec;
427 // mean energy
428 double emean = h->GetMean();
429 KVDetector* det = gMultiDetArray->GetDetector(name);
430 double fluence = rate / det->GetEntranceWindowSurfaceArea();
431 double dissipation = emean * rate / det->GetEntranceWindowSurfaceArea();
432 count_rates.push_back(
433 count_rate(name, rate, emean, map->GetMean(), map->GetMean(2), fluence, dissipation, intXsec)
434 );
435 fRates[name.Data()] = KVElasticCountRate(rate, emean, intXsec, fluence, dissipation);
436 }
437 }
438 std::sort(count_rates.begin(), count_rates.end(), compare_count_rates);
439
440 for (std::vector<count_rate>::iterator it = count_rates.begin(); it != count_rates.end(); ++it) {
441 it->print();
442 }
443}
444
445
446
449
451{
452 // Print mean energy deposit & counting rate for given beam intensity in particles per second
454 nl.SetName("Generated by KVElasticCountRates::PutResultsInList method");
455 KVDatime dt;
456 nl.SetTitle(dt.AsSQLString());
457 TIter it(&fHistos);
458 TH1F* h;
459 fRates.clear();
460 nl.SetValue("Intensity(pps)", beam_intensity);
461
462 std::vector<count_rate> count_rates;
463
464 while ((h = (TH1F*)it())) {
465 TString name = h->GetName();
466 if (!name.EndsWith("_dW") && !name.EndsWith("_map")) {
467 TH2F* map = (TH2F*)fHistos.FindObject(name + "_map");
468 // integrated cross-section
469 double intXsec = h->Integral() * fVolume / fNtirages;
470 // counting rate
471 double rate = fAtomicDensity * beam_intensity * intXsec;
472 // mean energy
473 double emean = h->GetMean();
474 KVDetector* det = gMultiDetArray->GetDetector(name);
475 double fluence = rate / det->GetEntranceWindowSurfaceArea();
476 double dissipation = emean * rate / det->GetEntranceWindowSurfaceArea();
477 count_rates.push_back(
478 count_rate(name, rate, emean, map->GetMean(), map->GetMean(2), fluence, dissipation, intXsec)
479 );
480 fRates[name.Data()] = KVElasticCountRate(rate, emean, intXsec, fluence, dissipation);
481 }
482 }
483 std::sort(count_rates.begin(), count_rates.end(), compare_count_rates);
484
485 for (std::vector<count_rate>::iterator it = count_rates.begin(); it != count_rates.end(); ++it) {
486 it->PutInList(nl);
487 }
488 return nl;
489}
490
491
492//__________________________________________________________________//
493
int Int_t
#define d(i)
#define f(i)
#define c(i)
#define e(i)
char Char_t
constexpr Bool_t kFALSE
double Double_t
#define N
winID h TVirtualViewer3D TVirtualGLPainter p
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 b
char name[80]
char * Form(const char *fmt,...)
Relativistic binary kinematics calculator.
Definition KV2Body.h:166
void SetTarget(const KVNucleus &)
Set target for reaction.
Definition KV2Body.cpp:314
void SetEDiss(Double_t ex)
Definition KV2Body.h:246
void SetProjectile(const KVNucleus &)
Set projectile for reaction.
Definition KV2Body.cpp:339
Double_t GetXSecRuthLab(Double_t ThetaLab_Proj, Int_t OfNucleus=3) const
Definition KV2Body.cpp:1279
void SetOutgoing(const KVNucleus &proj_out)
Definition KV2Body.cpp:397
Double_t GetELab(Double_t ThetaCM, Int_t OfNucleus) const
Definition KV2Body.h:291
void CalculateKinematics()
Definition KV2Body.cpp:677
Extension of TDatime to handle various useful date formats.
Definition KVDatime.h:33
Base class for detector geometry description.
Definition KVDetector.h:160
virtual void Clear(Option_t *opt="")
virtual Double_t GetEntranceWindowSurfaceArea()
Return surface area of first layer of detector in cm2.
Calculate elastic scattering count rates in multidetector arrays.
Int_t fIntLayer
index of interaction layer in multilayer target
void PrintResults(Double_t beam_intensity=1.e+07)
Print mean energy deposit & counting rate for given beam intensity in particles per second.
Double_t fEnergy
energy of projectile
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 *)
KVNucleus fProj
scattered nucleus
TH2F * fGlobalMap
global map of cross-section vs. theta vs. phi
void SetProjectile(const Char_t *nuc, Double_t e_sur_a)
Set projectile and beam energy [MeV/nucleon] & direction.
Int_t fBinE
Number of bins of the Energy histogram.
void SetEbinning(Int_t nbins=500)
KVHashList fHistos
histograms for all hit detectors
KV2Body * fKinematics
kinematics calculation
void CalculateScattering(Int_t N=10000)
Double_t fVolume
volume factor for MC integration
void SetRun(Int_t run)
Set detector parameters, target, etc. for run.
TH1F * fDepth
depth of scattering point in target
std::map< std::string, KVElasticCountRate > fRates
Bool_t fMultiLayer
kTRUE for multilayer target
void SetTargetScatteringLayer(const Char_t *name)
virtual ~KVElasticCountRates()
Destructor.
Double_t fExx
excited state of target nucleus
KVPosition fAngularRange
angular range in which to scatter
TH1F * fTheta
angle of scattered particle
TVector3 fBeamDirection
beam direction vector
KVTarget * fTarget
target for current run
Double_t fAtomicDensity
number of atoms per barn (10^-24 cm2) in target
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:94
Double_t GetZ() const
Double_t GetMass() const
virtual KVNameValueList * DetectParticle(KVNucleus *part)
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)
KVTarget * GetTarget()
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)
virtual void Clear(Option_t *opt="")
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:126
Bool_t IsDefined() const
Definition KVNucleus.h:202
void SetZandA(Int_t z, Int_t a)
Set atomic number and mass number.
void SetTheta(Double_t theta)
Definition KVParticle.h:693
void SetPhi(Double_t phi)
Definition KVParticle.h:697
Double_t GetE() const
Definition KVParticle.h:652
void SetMomentum(const TVector3 *v)
Definition KVParticle.h:542
KVNameValueList * GetParameters() const
Definition KVParticle.h:815
void SetEnergy(Double_t e)
Definition KVParticle.h:599
virtual void GetRandomAngles(Double_t &th, Double_t &ph, Option_t *t="isotropic")
Double_t GetThetaMax() const
Definition KVPosition.h:154
virtual TObject * FindObject(const char *name) const
virtual void Add(TObject *obj)
virtual void Delete(Option_t *option="")
Extension of ROOT TString class which allows backwards compatibility with ROOT v3....
Definition KVString.h:73
void SetIncoming(Bool_t r=kTRUE)
Definition KVTarget.h:217
virtual void DetectParticle(KVNucleus *, TVector3 *norm=0)
Definition KVTarget.cpp:485
Int_t NumberOfLayers() const
Definition KVTarget.h:166
TVector3 & GetInteractionPoint(KVParticle *part=0)
Definition KVTarget.cpp:767
void SetRandomized(Bool_t r=kTRUE)
Definition KVTarget.h:250
KVMaterial * GetLayer(TVector3 &depth)
Definition KVTarget.cpp:272
Double_t GetTotalEffectiveThickness(KVParticle *part=0)
Definition KVTarget.cpp:412
void SetInteractionLayer(Int_t ilayer, TVector3 &dir)
Definition KVTarget.cpp:809
void SetOutgoing(Bool_t r=kTRUE)
Definition KVTarget.h:237
Double_t GetAtomsPerCM2() const
virtual UInt_t GetUnits() const;
Definition KVTarget.cpp:926
Int_t GetLayerIndex(TVector3 &depth)
Definition KVTarget.cpp:292
const char * AsSQLString() const
virtual Double_t GetMean(Int_t axis=1) 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)
Double_t Angle(const TVector3 &v) const
virtual void SetTitle(const char *title="")
const char * GetName() const override
virtual void SetName(const char *name)
const char * Data() const
Double_t z() const
RVec< PromoteType< T > > cos(const RVec< T > &v)
RVec< PromoteType< T > > sin(const RVec< T > &v)
const Int_t n
TH1 * h
void Warning(const char *location, const char *fmt,...)
constexpr Double_t DegToRad()
Double_t Abs(Double_t d)
Double_t Max(Double_t a, Double_t b)
constexpr Double_t RadToDeg()
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)
TArc a
ClassImp(TPyArg)