BmnRoot
Loading...
Searching...
No Matches
BmnParticleEqualizer.h
Go to the documentation of this file.
1// -------------------------------------------------------------------------
2// ----- BmnParticleEqualizer header file -----
3// ----- Created 01/05/2024 by K. Gertsenberger -----
4// -------------------------------------------------------------------------
5
12#ifndef BMNPARTICLEEQUALIZER_H
13#define BMNPARTICLEEQUALIZER_H
14
15#include "FairGenerator.h"
16#include "FairPrimaryGenerator.h"
17#include "TF1.h"
18#include "TH1D.h"
19
20#include <map>
21#include <vector>
22using namespace std;
23
24class BmnParticleEqualizer : public FairGenerator
25{
26 public:
28 BmnParticleEqualizer(TString hist_file_name);
35 BmnParticleEqualizer(TString hist_file_name,
36 vector<Int_t> pdg_codes,
37 Double_t min_p = 0.0,
38 Double_t max_p = 4.0,
39 Int_t intervals = 20);
40
43
45 Bool_t Init();
49 virtual Bool_t ReadEvent(FairPrimaryGenerator* primGen);
50
52 {
53 // Default 'stParticleInfo' constructor for storing particle distrubutons
54 stParticleInfo(int pdg_key, Double_t min_p, Double_t max_p, Int_t interval_count);
55 // 'stParticleInfo' destructor, empty becausse of using unique pointers
57
58 vector<Int_t> fParticleNumber;
59
60 unique_ptr<TH1D> fPHist;
61 unique_ptr<TH1D> fPtHist;
62 unique_ptr<TH1D> fYHist;
63 unique_ptr<TH1D> fEtaHist;
64 unique_ptr<TH1D> fThetaHist;
65 unique_ptr<TH1D> fPhiHist;
66
67 vector<unique_ptr<TH1D>>
69 vector<unique_ptr<TH1D>>
71
72 Double_t fY0;
73 Double_t fSigma;
74 Double_t fT;
75 };
76
77 // static function to form sample histograms for equalizing the particle distributions
78 static int ProduceSampleHistograms(TString input_list_file,
79 TString output_histo_file = "$VMCWORKDIR/macro/recotools/particle_hists.root");
80
81 // fit P. Pt, Y, Eta histograms to generate the values according to the analytical functions
82 static int FitPHistogram(unique_ptr<TH1D>& hP, int pdg_key, Double_t p_min, Double_t p_max);
83 static int FitPtHistogram(unique_ptr<TH1D>& hPt,
84 int pdg_key,
85 Double_t pdg_mass,
86 Double_t pt_min,
87 Double_t pt_max,
88 Double_t& t_output);
89 static int FitYHistogram(unique_ptr<TH1D>& hY, int pdg_key, Double_t& y0_output, Double_t& sigma_putput);
90 static int FitEtaHistogram(unique_ptr<TH1D>& hEta, int pdg_key, Double_t eta_min, Double_t eta_max);
91
93 // static function to fit the obtained sample histograms (to show whether the results are appropriate)
94 static void FitSampleHistograms(TString input_histo_file = "$VMCWORKDIR/macro/recotools/particle_hists.root");
95 // static function to show the obtained sample histograms with obtained fit (to check whether the results are
96 // appropriate)
97 static int ShowSampleHistograms(TString input_histo_file = "$VMCWORKDIR/macro/recotools/particle_hists.root");
98 // draw histogram with multiplicity of different particle types for all the momentum intervals
99 // as well as P, Pt, Y, Eta, Theta, Phi distributions after the equalizing
100 // it work both for SIM and DST files listed in the given text file ('input_list_file' var)
101 static int ShowResultDistributions(TString input_list_file);
102
103 private:
104 static Double_t fPMin;
105 static Double_t fPMax;
106 static Int_t fIntervalCount;
107 Double_t fIntervalStep;
108
109 // Whether correct input data for the generator is obtained (if false, it will be skipped)
110 Bool_t isActive;
111 // Particle PDG codes to be equalized
112 static vector<Int_t> fParticles;
113 // map with the full information on all the particle types
114 map<Int_t, unique_ptr<stParticleInfo>> fParticleInfos;
115 // Maximum number of a particle type among all particles (divided by P intervals -> vector)
116 vector<int> fMaxParticles;
117
118 // Get all sample histograms for all the particle types from the 'input_histo_file' file
119 // and calculate fit parameters (if necessary) according to the given histograms
120 Bool_t ReadSampleHistograms(TString input_histo_file = "$VMCWORKDIR/macro/recotools/particle_hists.root");
121
122 // define PDG code of a particle by the given mass
123 // among particle types set by 'vecParticlePDG' vector or in the TDatabasePDG if the vector is empty
124 // tolerance: relative if tolerance > 0, absolute (in GeV) in case of tolerance <= 0
125 // return PDG code or 0 if no corresponding particle was found
126 static Int_t GetPDGCodeByMass(Double_t mass, Double_t tolerance = 1.e-6, const vector<int>& vecParticlePDG = {});
127
128 ClassDef(BmnParticleEqualizer, 1);
129};
130
131#endif
static int FitEtaHistogram(unique_ptr< TH1D > &hEta, int pdg_key, Double_t eta_min, Double_t eta_max)
static int ShowSampleHistograms(TString input_histo_file="$VMCWORKDIR/macro/recotools/particle_hists.root")
static void FitSampleHistograms(TString input_histo_file="$VMCWORKDIR/macro/recotools/particle_hists.root")
static int FitPHistogram(unique_ptr< TH1D > &hP, int pdg_key, Double_t p_min, Double_t p_max)
static int FitPtHistogram(unique_ptr< TH1D > &hPt, int pdg_key, Double_t pdg_mass, Double_t pt_min, Double_t pt_max, Double_t &t_output)
virtual Bool_t ReadEvent(FairPrimaryGenerator *primGen)
static int ShowResultDistributions(TString input_list_file)
static int ProduceSampleHistograms(TString input_list_file, TString output_histo_file="$VMCWORKDIR/macro/recotools/particle_hists.root")
static int FitYHistogram(unique_ptr< TH1D > &hY, int pdg_key, Double_t &y0_output, Double_t &sigma_putput)
STL namespace.
vector< unique_ptr< TH1D > > fEtaHistVector
Eta histograms for the particle type (divided by momentum intervals -> vector)
unique_ptr< TH1D > fPHist
PDG -> P histogram (TH1D) for the particle type.
Double_t fT
Temperature in the P distribution.
vector< unique_ptr< TH1D > > fPHistVector
P histograms for the particle type (divided by momentum intervals -> vector)
unique_ptr< TH1D > fPhiHist
PDG -> Phi (azimuth angle) histogram (TH1D) for the particle type.
unique_ptr< TH1D > fEtaHist
PDG -> Eta (pseudorapidity) histogram (TH1D) for the particle type.
unique_ptr< TH1D > fPtHist
PDG -> Pt histogram (TH1D) for the particle type.
unique_ptr< TH1D > fYHist
PDG -> Y (rapidity) histogram (TH1D) for the particle type.
Double_t fSigma
Sigma in the rapidity distribution.
unique_ptr< TH1D > fThetaHist
PDG -> Theta (polar angle) histogram (TH1D) for the particle type.
vector< Int_t > fParticleNumber
PDG -> Particle Count (divided by momentum intervals -> vector)