KaliVeda
Toolkit for HIC analysis
KVElasticCountRates Class Reference

Detailed Description

Calculate elastic scattering count rates in the detectors of multidetector arrays ,.

Created by KVClassFactory on Fri Nov 20 2015 Author: John Frankland

Use this class to calculate the count rates & energy losses of elastically scattered nuclei in the detectors of a multidetector array corresponding to a given angular range.

For details on setting up and fine-tuning the simulation, see class KVElasticScatter.

Running the calculation, obtaining the results

To perform a calculation, the minimum you need to do is:

KVMultiDetArray::MakeMultiDetector("name_of_dataset"); // build geometry
KVElasticCountRates ecr(run_number); // choose run => define projectile, target etc.
ecr.SetAngularRange(theta_min,theta_max,phi_min,phi_max); // all detectors in range will be sampled
ecr.CaclulateScattering(10000);
Calculate elastic scattering count rates in the detectors of multidetector arrays ,...
static KVMultiDetArray * MakeMultiDetector(const Char_t *dataset_name, Int_t run=-1, TString classname="KVMultiDetArray", KVExpDB *db=nullptr)

Make sure the number of sampled points (default=10000) is large enough for accurate determination of count rates.

At the end of calculation, we print infos for all detectors hit by elastic particles, these are the count rates for a nominal beam intensity of 10**7 particles per second, the mean energy loss in the detector, the total integrated cross-section per detector in barn, the number of particles per second per unit area hitting each detector, and the total energy per second per unit area deposited in each detector:

esa.CalculateScattering(10000);
SI_0101 : N= 204.42/sec. <E>= 890.3 MeV Tot.Xsec=6.687E+01 barn fluence=2.711E+01/sec./cm**2 dissip.=2.413E+04 MeV/sec./cm**2
SI_0201 : N= 61.37/sec. <E>= 709.7 MeV Tot.Xsec=2.008E+01 barn fluence=1.413E+01/sec./cm**2 dissip.=1.002E+04 MeV/sec./cm**2
SI_0301 : N= 15.29/sec. <E>= 706.4 MeV Tot.Xsec=5.002E+00 barn fluence=2.390E+00/sec./cm**2 dissip.=1.688E+03 MeV/sec./cm**2
SI_0401 : N= 6.33/sec. <E>= 699.2 MeV Tot.Xsec=2.070E+00 barn fluence=1.582E+00/sec./cm**2 dissip.=1.106E+03 MeV/sec./cm**2
#define N
constexpr Double_t E()
const long double cm
Definition: KVUnits.h:69
const long double MeV
energies
Definition: KVUnits.h:91

If you want to see the count rates for a different beam intensity, call

ecr.PrintResults(5.e+06)

with the required number of projectiles per second.

For each detector hit we fill three histograms:

OBJ: TH1F SI_0223 Eloss in SI_0223 : 0 at: 0x10898980
OBJ: TH1F SI_0223_dW Solid angle of SI_0223 : 0 at: 0x10899540
OBJ: TH2F SI_0223_map Map of SI_0223 : 0 at: 0x10e2e7d0
Option_t Option_t TPoint TPoint angle
auto Map(Args &&... args)

Each histogram can be obtained using

ecr.GetHistos().FindObject("SI_0223_map")

The different informations for each detector can be obtained using

ecr.GetDetector("SI_0201").count_rate
(double)6.13702995974404857e+01
ecr.GetDetector("SI_0201").mean_energy
(double)7.09659910827807721e+02
ecr.GetDetector("SI_0201").intXsec
(double)2.00758282694385635e+01

See KVElasticCountRate class for details.

Warning
The values returned correspond to the last beam intensity value given to PrintResults().

Definition at line 115 of file KVElasticCountRates.h.

#include <KVElasticCountRates.h>

Inheritance diagram for KVElasticCountRates:

Public Member Functions

 KVElasticCountRates (int run)
 
void FillHistograms (const KVNameValueList *, double theta, double phi, double xsec)
 
KVElasticCountRate GetDetector (const std::string &name)
 
const KVHashListGetHistos () const
 
void PrintResults (Double_t beam_intensity=1.e+07)
 Print mean energy deposit & counting rate for given beam intensity in particles per second. More...
 
KVNameValueList PutResultsInList (Double_t beam_intensity=1.e+07)
 Print mean energy deposit & counting rate for given beam intensity in particles per second. More...
 
void SetAngularRange (Double_t theta_min=0, Double_t theta_max=180, Double_t phi_min=0, Double_t phi_max=360)
 
void SetBeamIntensity (Double_t bi)
 
- Public Member Functions inherited from KVElasticScatter
 KVElasticScatter (Int_t run)
 Default constructor. More...
 
virtual ~ KVElasticScatter ()
 
void CalculateScattering (Int_t N)
 Perform scattering 'N' times. More...
 
TH1FGetDepth ()
 Return pointer to histogram of 'depth' of scattering point in target (in mg/cm2) More...
 
Int_t GetEbinning (void)
 Returns the number of bins of the GetEnergy histogram. More...
 
TH1FGetEnergy ()
 Return pointer to energy loss histogram for chosen detector (in MeV) More...
 
TH1FGetEnergy (const Char_t *label)
 Energy loss in detector with given 'label' through which scattered particle passes. More...
 
TH1FGetEnergy (Int_t index)
 
auto GetKinematics ()
 
Int_t GetNDets () const
 Returns the number of detectors crossed by the scattered particle. More...
 
void GetTargetNucleusFromTargetLayer (Bool_t yes=kTRUE)
 
TH1FGetTheta ()
 Return pointer to polar angle distribution of scattered particle (in degrees) More...
 
void Print () const
 Print details of calculation to be performed. More...
 
void SetBeamDirection (const TVector3 &beam_dir)
 
void SetDetector (const Char_t *det)
 
void SetEbinning (Int_t nbins)
 
void SetExcitation (Double_t ex)
 
void SetKinematicSolution (KV2Body::kinematic_solution ik)
 
void SetNucleusOfInterest (KV2Body::nucleus_of_interest noi)
 
void SetOutGoingProjectileLike (const KVNucleus &opl)
 
void SetTargetNucleusForScattering (const KVNucleus &tn)
 
void SetTargetScatteringLayer (const Char_t *name)
 

Private Member Functions

void detect_particle_fill_histograms (KVNucleus *ejectile, double theta, double phi, double xsec) override
 now detect particle in array More...
 
void end_of_run () override
 
std::pair< double, double > get_random_angles_for_scattering (const KV2Body &) override
 set random direction of outgoing projectile More...
 
bool initial_checks_and_reset () override
 
void reset_before_new_scattering () override
 

Private Attributes

KVDetectionSimulator det_sim
 
KVPosition fAngularRange
 angular range in which to scatter More...
 
KVHashList fArrayHistos
 histograms for all hit detectors More...
 
Double_t fBeamIntensity = 1.e+7
 
TH2FfGlobalMap
 global map of cross-section vs. theta vs. phi More...
 
std::map< std::string, KVElasticCountRatefRates
 
Double_t fVolume = 0
 
Double_t phi
 
Double_t theta
 

Constructor & Destructor Documentation

◆ KVElasticCountRates()

KVElasticCountRates::KVElasticCountRates ( int  run)
inline

Definition at line 132 of file KVElasticCountRates.h.

Member Function Documentation

◆ detect_particle_fill_histograms()

void KVElasticCountRates::detect_particle_fill_histograms ( KVNucleus ejectile,
double  theta,
double  phi,
double  xsec 
)
overrideprivatevirtual

now detect particle in array

Reimplemented from KVElasticScatter.

Definition at line 80 of file KVElasticCountRates.cpp.

◆ end_of_run()

void KVElasticCountRates::end_of_run ( )
overrideprivatevirtual

Reimplemented from KVElasticScatter.

Definition at line 91 of file KVElasticCountRates.cpp.

◆ FillHistograms()

void KVElasticCountRates::FillHistograms ( const KVNameValueList dets,
double  theta,
double  phi,
double  xsec 
)

parse the list dets fill histograms with energy loss for all detectors clear the detector energy losses delete the list

Definition at line 104 of file KVElasticCountRates.cpp.

◆ get_random_angles_for_scattering()

std::pair< double, double > KVElasticCountRates::get_random_angles_for_scattering ( const KV2Body )
overrideprivatevirtual

set random direction of outgoing projectile

Reimplemented from KVElasticScatter.

Definition at line 67 of file KVElasticCountRates.cpp.

◆ GetDetector()

KVElasticCountRate KVElasticCountRates::GetDetector ( const std::string &  name)
inline

Definition at line 150 of file KVElasticCountRates.h.

◆ GetHistos()

const KVHashList& KVElasticCountRates::GetHistos ( ) const
inline

Definition at line 142 of file KVElasticCountRates.h.

◆ initial_checks_and_reset()

bool KVElasticCountRates::initial_checks_and_reset ( )
overrideprivatevirtual

Reimplemented from KVElasticScatter.

Definition at line 43 of file KVElasticCountRates.cpp.

◆ PrintResults()

void KVElasticCountRates::PrintResults ( Double_t  beam_intensity = 1.e+07)

Print mean energy deposit & counting rate for given beam intensity in particles per second.

Definition at line 202 of file KVElasticCountRates.cpp.

◆ PutResultsInList()

KVNameValueList KVElasticCountRates::PutResultsInList ( Double_t  beam_intensity = 1.e+07)

Print mean energy deposit & counting rate for given beam intensity in particles per second.

Definition at line 243 of file KVElasticCountRates.cpp.

◆ reset_before_new_scattering()

void KVElasticCountRates::reset_before_new_scattering ( )
overrideprivatevirtual

Reimplemented from KVElasticScatter.

Definition at line 57 of file KVElasticCountRates.cpp.

◆ SetAngularRange()

void KVElasticCountRates::SetAngularRange ( Double_t  theta_min = 0,
Double_t  theta_max = 180,
Double_t  phi_min = 0,
Double_t  phi_max = 360 
)

Definition at line 33 of file KVElasticCountRates.cpp.

◆ SetBeamIntensity()

void KVElasticCountRates::SetBeamIntensity ( Double_t  bi)
inline

Definition at line 137 of file KVElasticCountRates.h.

Member Data Documentation

◆ det_sim

KVDetectionSimulator KVElasticCountRates::det_sim
private

Definition at line 124 of file KVElasticCountRates.h.

◆ fAngularRange

KVPosition KVElasticCountRates::fAngularRange
private

angular range in which to scatter

Definition at line 117 of file KVElasticCountRates.h.

◆ fArrayHistos

KVHashList KVElasticCountRates::fArrayHistos
private

histograms for all hit detectors

Definition at line 120 of file KVElasticCountRates.h.

◆ fBeamIntensity

Double_t KVElasticCountRates::fBeamIntensity = 1.e+7
private

Definition at line 122 of file KVElasticCountRates.h.

◆ fGlobalMap

TH2F* KVElasticCountRates::fGlobalMap
private

global map of cross-section vs. theta vs. phi

Definition at line 116 of file KVElasticCountRates.h.

◆ fRates

std::map<std::string, KVElasticCountRate> KVElasticCountRates::fRates
private

Definition at line 121 of file KVElasticCountRates.h.

◆ fVolume

Double_t KVElasticCountRates::fVolume = 0
private

Definition at line 118 of file KVElasticCountRates.h.

◆ phi

Double_t KVElasticCountRates::phi
private

Definition at line 119 of file KVElasticCountRates.h.

◆ theta

Double_t KVElasticCountRates::theta
private

Definition at line 119 of file KVElasticCountRates.h.