KaliVeda
Toolkit for HIC analysis
KVElasticScatter.h
1 /*
2 $Id: KVElasticScatter.h,v 1.6 2007/04/04 10:39:17 ebonnet Exp $
3 $Revision: 1.6 $
4 $Date: 2007/04/04 10:39:17 $
5 */
6 
9 
10 #ifndef __KVELASTICSCATTER_H
11 #define __KVELASTICSCATTER_H
12 
13 #include "TVector3.h"
14 #include "KVNameValueList.h"
15 #include "KVNucleus.h"
16 #include "KV2Body.h"
17 #include <optional>
18 
19 class TH1F;
20 class TList;
21 class KVGroup;
22 class KVDetector;
23 class KVTarget;
24 class TObjArray;
25 
191 
192  TH1F* fDepth = nullptr; //depth of scattering point in target
193  TH1F* fTheta = nullptr; //angle of scattered particle
194  Int_t fBinE = 100; //Number of bins of the Energy histogram
195  KV2Body* fKinematics = nullptr; //kinematics calculation
196  KVDetector* fDetector = nullptr; //detector where particle will be detected
199  KVTarget* fTarget = nullptr; //target for current run
200  Bool_t fMultiLayer = false; //kTRUE for multilayer target
201  Bool_t fGetTargFromLay = false; //kTRUE if target nucleus can change from layer to layer
202  TVector3 fBeamDirection{0, 0, 1}; //beam direction vector
203  Int_t fIntLayer = 0; //index of interaction layer in multilayer target
206  Double_t fExx = 0.0; // excitation of target
207  std::optional<KVNucleus> fOutgoingPL;
208  std::optional<KVNucleus> fTargetNuc;
209  void SetRun(Int_t run);
210  KV2Body::nucleus_of_interest fNucleusOfInterest = KV2Body::nucleus_of_interest::projectile_like;
211  KV2Body::kinematic_solution fKineSol = KV2Body::kinematic_solution::high_E_branch;
212  Double_t fFWHM = 0.; // detector resolution at some reference energy
213  Double_t fFWHM_ref = 0; // the reference energy
214 
215  double get_det_resolution(double E)
216  {
217  return fFWHM*TMath::Sqrt(E*fFWHM_ref);
218  }
219 
220 protected:
223  virtual bool initial_checks_and_reset();
224  virtual void reset_before_new_scattering() {}
225  virtual std::pair<double, double> get_random_angles_for_scattering(const KV2Body& scattering_kinematics);
226  virtual void detect_particle_fill_histograms(KVNucleus* ejectile, double theta, double phi, double xsec);
227  virtual void end_of_run() {}
228 
229 public:
230  KVElasticScatter(Int_t run);
232 
233  void SetBeamDirection(const TVector3& beam_dir)
234  {
236  fBeamDirection = beam_dir;
237  }
238  void SetTargetScatteringLayer(const Char_t* name);
240  {
249  fNucleusOfInterest = noi;
250  }
252  {
261  fKineSol = ik;
262  }
263 
264  void SetDetectorResolution(double fwhm, double eref)
265  {
269  fFWHM = fwhm;
270  fFWHM_ref = eref;
271  }
272  void SetEbinning(Int_t nbins);
273  void SetDetector(const Char_t* det);
274  void CalculateScattering(Int_t N);
276  {
282  fExx = ex;
283  }
285  {
286  fOutgoingPL = opl;
287  fKinematics->SetOutgoing(opl);
289  }
291  {
297  fGetTargFromLay = yes;
298  }
300  {
302  fTargetNuc = tn;
303  fKinematics->SetTarget(tn);
305  }
307  {
311  }
312 
315  {
316  return fDepth;
317  }
320  {
321  return GetEnergy(fNDets - 1);
322  }
325  {
326  return fTheta;
327  }
328  TH1F* GetEnergy(const Char_t* label);
329  TH1F* GetEnergy(Int_t index);
331  Int_t GetNDets() const
332  {
333  return fNDets;
334  }
337  {
338  return fBinE;
339  }
341  {
342  return fKinematics;
343  }
344 
345  void Print() const;
346 
347  ClassDef(KVElasticScatter, 0) //Calculate elastic scattering spectra in multidetector arrays
348 };
349 
350 #endif
int Int_t
bool Bool_t
char Char_t
double Double_t
#define ClassDef(name, id)
Relativistic binary kinematics calculator.
Definition: KV2Body.h:178
void SetTarget(const KVNucleus &)
Set target for reaction.
Definition: KV2Body.cpp:328
void SetProjectile(const KVNucleus &)
Set projectile for reaction.
Definition: KV2Body.cpp:353
void SetOutgoing(const KVNucleus &proj_out)
Definition: KV2Body.cpp:411
nucleus_of_interest
Definition: KV2Body.h:263
void CalculateKinematics()
Definition: KV2Body.cpp:693
kinematic_solution
Definition: KV2Body.h:269
Base class for detector geometry description, interface to energy-loss calculations.
Definition: KVDetector.h:160
Calculate elastic scattering spectra in specific detectors of a multidetector array ,...
virtual std::pair< double, double > get_random_angles_for_scattering(const KV2Body &scattering_kinematics)
void SetNucleusOfInterest(KV2Body::nucleus_of_interest noi)
void SetDetector(const Char_t *det)
KV2Body::nucleus_of_interest fNucleusOfInterest
Int_t fNDets
number of aligned detectors
virtual void end_of_run()
std::optional< KVNucleus > fOutgoingPL
optional different nucleus for outgoing projectile-like
virtual ~ KVElasticScatter()
void SetTargetNucleusForScattering(const KVNucleus &tn)
virtual bool initial_checks_and_reset()
void SetProjectileNucleus(const KVNucleus &pn)
void SetBeamDirection(const TVector3 &beam_dir)
virtual void reset_before_new_scattering()
std::optional< KVNucleus > fTargetNuc
optional different target nucleus for scattering
KVDetector * fDetector
TH1F * GetEnergy()
Return pointer to energy loss histogram for chosen detector (in MeV)
void CalculateScattering(Int_t N)
Perform scattering 'N' times.
void Print() const
Print details of calculation to be performed.
Int_t GetNDets() const
Returns the number of detectors crossed by the scattered particle.
void SetDetectorResolution(double fwhm, double eref)
void SetRun(Int_t run)
TList fAlignedDetectors
all aligned detectors
void SetKinematicSolution(KV2Body::kinematic_solution ik)
TObjArray fHistos
energy loss histograms for all hit detectors
KVElasticScatter(Int_t run)
Default constructor.
void SetTargetScatteringLayer(const Char_t *name)
Double_t fAtomicDensity
number of atoms per barn (10^-24 cm2) in target
Int_t GetEbinning(void)
Returns the number of bins of the GetEnergy histogram.
KVNameValueList fDetInd
detector type-index association
virtual void detect_particle_fill_histograms(KVNucleus *ejectile, double theta, double phi, double xsec)
KV2Body::kinematic_solution fKineSol
double get_det_resolution(double E)
TH1F * GetDepth()
Return pointer to histogram of 'depth' of scattering point in target (in mg/cm2)
void SetOutGoingProjectileLike(const KVNucleus &opl)
void SetEbinning(Int_t nbins)
TH1F * GetTheta()
Return pointer to polar angle distribution of scattered particle (in degrees)
void SetExcitation(Double_t ex)
void GetTargetNucleusFromTargetLayer(Bool_t yes=kTRUE)
Group of detectors which can be treated independently of all others in array.
Definition: KVGroup.h:20
Handles lists of named parameters with different types, a list of KVNamedParameter objects.
Description of properties and kinematics of atomic nuclei.
Definition: KVNucleus.h:123
Calculation/correction of energy losses of particles through an experimental target.
Definition: KVTarget.h:128
Double_t ex[n]
constexpr Double_t E()
Double_t Sqrt(Double_t x)