KaliVeda
Toolkit for HIC analysis
Loading...
Searching...
No Matches
KVElasticCountRates.h
1
3
4#ifndef __KVElasticCountRates_H
5#define __KVElasticCountRates_H
6
7#include "TVector3.h"
8#include "KVNameValueList.h"
9#include "KVPosition.h"
10#include "TH2.h"
11#include "KVNucleus.h"
12#include "KV2Body.h"
13#include "KVTarget.h"
14#include <map>
15
140 double intXsec; // integrated X-section (barn)
141 double fluence; // ions/sec./cm**2
142 double dissipation; // MeV/sec./cm**2
143 KVElasticCountRate(double c = 0, double e = 0, double i = 0, double f = 0, double d = 0)
146
147 ClassDef(KVElasticCountRate, 0) //Elastic scattering rate information for detector
148};
149
151
156
159
161
163
171
174
176
178
179 std::map<std::string, KVElasticCountRate> fRates;
180
181public:
182
183 KVElasticCountRates(Double_t theta_min = 0, Double_t theta_max = 180, Double_t phi_min = 0, Double_t phi_max = 360, const TVector3& beam_dir = {0, 0, 1});
184 virtual ~ KVElasticCountRates();
185
186 void SetRun(Int_t run);
187 void SetProjectile(const Char_t* nuc, Double_t e_sur_a);
188 void SetTargetScatteringLayer(const Char_t* name);
189 void SetEbinning(Int_t nbins = 500);
190
191 void CalculateScattering(Int_t N = 10000);
193 {
198 fExx = ex;
199 }
200
202 {
204
205 return fDepth;
206 }
208 {
210
211 return fTheta;
212 }
214 {
216
217 return fBinE;
218 }
219
221 const KVHashList& GetHistos() const
222 {
223 return fHistos;
224 };
225
226 void PrintResults(Double_t beam_intensity = 1.e+07);
227 KVNameValueList PutResultsInList(Double_t beam_intensity = 1.e+07);
228
229 KVElasticCountRate GetDetector(const std::string& name)
230 {
231 return fRates[name];
232 }
233
234 ClassDef(KVElasticCountRates, 1)//Calculate elastic scattering count rates in multidetector arrays
235};
236
237#endif
int Int_t
#define d(i)
#define f(i)
#define c(i)
#define e(i)
bool Bool_t
char Char_t
double Double_t
#define ClassDef(name, id)
char name[80]
Relativistic binary kinematics calculator.
Definition KV2Body.h:166
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
KVElasticCountRate GetDetector(const std::string &name)
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)
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
void SetTargetExcitedState(Double_t ex)
Double_t fAtomicDensity
number of atoms per barn (10^-24 cm2) in target
const KVHashList & GetHistos() const
Extended version of ROOT THashList.
Definition KVHashList.h:29
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Description of properties and kinematics of atomic nuclei.
Definition KVNucleus.h:126
Base class used for handling geometry in a multidetector array.
Definition KVPosition.h:91
Calculation/correction of energy losses of particles through an experimental target.
Definition KVTarget.h:127
Double_t ex[n]
Utility class used by KVElasticCountRates to store results.
KVElasticCountRate(double c=0, double e=0, double i=0, double f=0, double d=0)
Utility class used by KVElasticCountRates.