StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
dEdxParameterization.h
1 
2 #ifndef dEdxParameterization_h
3 #define dEdxParameterization_h
4 
5 #include "TROOT.h"
6 #include "TChain.h"
7 #include "TFile.h"
8 #include "TH1.h"
9 #include "TH2.h"
10 #include "TH3.h"
11 #include "TProfile2D.h"
12 #include "TString.h"
13 #include "StPidParticleDefinition.h"
15  private:
16  TString fTag;
17  TProfile2D *fP;
18  TProfile2D *fA;
19  TProfile2D *fI70;
20  TProfile2D *fI60;
21  TProfile2D *fD;
22  TProfile2D *fRms;
23  TProfile2D *fW;
24  TH3D *fPhi;
25  Int_t fnBins[3];
26  Double_t fbinW[3];
27  TAxis *fAXYZ[3];
28  Double_t fMostProbableZShift;
29  Double_t fAverageZShift;
30  Double_t fI70Shift;
31  Double_t fI60Shift;
32  Double_t fbgL10min;
33  Double_t fbgL10max;
34  Double_t fdxL2min;
35  Double_t fdxL2max;
36  Double_t fzmin;
37  Double_t fzmax;
38  TH1D *fTrs[KPidParticles+1][6];
39  public :
40  dEdxParameterization(const Char_t *Tag="p10", Int_t keep3D = 0,
41  const Double_t MostProbableZShift = 0,
42  const Double_t AverageZShift = 0,
43  const Double_t I70Shift = 1,
44  const Double_t I60Shift = 1);
45  virtual ~dEdxParameterization();
46  Double_t MostProbableZCorrection(Double_t log10bg);
47  Double_t I70Correction(Double_t log10bg);
48  Double_t GetMostProbableZ(Double_t log10bg, Double_t log2dx) {
49  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
50  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
51  return fMostProbableZShift+fP->Interpolate(log10bg,log2dx);
52  }
53  Double_t GetAverageZ(Double_t log10bg, Double_t log2dx) {
54  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
55  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
56  return fAverageZShift+MostProbableZCorrection(log10bg)+fA->Interpolate(log10bg,log2dx);
57  }
58  Double_t GetRmsZ(Double_t log10bg, Double_t log2dx) {
59  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
60  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
61  return fRms->Interpolate(log10bg,log2dx);
62  }
63  Double_t GetI70(Double_t log10bg, Double_t log2dx) {
64  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
65  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
66  return fI70Shift*fI70->Interpolate(log10bg,log2dx);
67  }
68  Double_t GetI60(Double_t log10bg, Double_t log2dx) {
69  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
70  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
71  return fI60Shift*fI60->Interpolate(log10bg,log2dx);
72  }
73  Double_t GetMostProbabledEdx(Double_t log10bg, Double_t log2dx) {
74  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
75  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
76  return fD->Interpolate(log10bg,log2dx);
77  }
78  Double_t GetdEdxWidth(Double_t log10bg, Double_t log2dx) {
79  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
80  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
81  return fW->Interpolate(log10bg,log2dx);
82  }
83  Double_t GetMostProbableZM(Double_t log10bg, Double_t log2dx) {
84  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
85  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
86  return MostProbableZCorrection(log10bg)+GetMostProbableZ(log10bg,log2dx);
87  }
88  Double_t GetAverageZM(Double_t log10bg, Double_t log2dx) {
89  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
90  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
91  return MostProbableZCorrection(log10bg)+GetAverageZ(log10bg,log2dx);
92  }
93  Double_t GetI70M(Double_t log10bg, Double_t log2dx) {
94  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
95  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
96  return I70Correction(log10bg)*GetI70(log10bg,log2dx);
97  }
98  Double_t GetProbability(Double_t log10bg, Double_t log2dx, Double_t z) {
99  log10bg = TMath::Max(fbgL10min, TMath::Min(fbgL10max,log10bg));
100  log2dx = TMath::Max(fdxL2min, TMath::Min(fdxL2max,log2dx));
101  z = TMath::Max(fzmin, TMath::Min(fzmax,z));
102  return fPhi->Interpolate(log10bg,log2dx,z);}
103  void Print();
104  const Char_t *Tag() const {return fTag.Data();}
105  const TProfile2D *P() const {return fP;}
106  const TProfile2D *A() const {return fA;}
107  const TProfile2D *I70() const {return fI70;}
108  const TProfile2D *I60() const {return fI60;}
109  const TProfile2D *D() const {return fD;}
110  const TProfile2D *Rms() const {return fRms;}
111  const TProfile2D *W() const {return fW;}
112  const TH3D *Phi() const {return fPhi;}
113  Double_t bgL10min() const {return fbgL10min;}
114  Double_t bgL10max() const {return fbgL10max;}
115  const TH1D *I70Trs( Int_t part = KPidParticles) const {return fTrs[part][0];} // Estimation for I70 from TpcRS
116  const TH1D *I70TrsB( Int_t part = KPidParticles) const {return fTrs[part][1];} // Estimation for I70 - Bichsel from TpcRS
117  const TH1D *I70TrsS( Int_t part = KPidParticles) const {return fTrs[part][2];} // Estimation for relative sigma beta*gamma dependence for I70 from TpcRS normalized to MIP
118  const TH1D *IfitTrs( Int_t part = KPidParticles) const {return fTrs[part][3];} // Estimation for Ifit from TpcRS
119  const TH1D *IfitTrsB(Int_t part = KPidParticles) const {return fTrs[part][4];} // Estimation for Ifit - Bichsel from TpcRS
120  const TH1D *IfitTrsS(Int_t part = KPidParticles) const {return fTrs[part][5];} // Estimation for relative sigma beta*gamma dependence for Ifit from TpcRS normalized to MIP
121  Double_t Get(const TH1D *hist, Double_t log10bg) const;
122  Double_t I70Trs (Int_t part, Double_t log10bg) const {return Get(fTrs[part][0], log10bg);} // Estimation for I70 from TpcRS
123  Double_t I70TrsB (Int_t part, Double_t log10bg) const {return Get(fTrs[part][1], log10bg);} // Estimation for I70 - Bichsel from TpcRS
124  Double_t I70TrsS (Int_t part, Double_t log10bg) const {return Get(fTrs[part][2], log10bg);} // Estimation for relative sigma beta*gamma dependence for I70 from TpcRS normalized to MIP
125  Double_t IfitTrs (Int_t part, Double_t log10bg) const {return Get(fTrs[part][3], log10bg);} // Estimation for Ifit from TpcRS
126  Double_t IfitTrsB(Int_t part, Double_t log10bg) const {return Get(fTrs[part][4], log10bg);} // Estimation for Ifit - Bichsel from TpcRS
127  Double_t IfitTrsS(Int_t part, Double_t log10bg) const {return Get(fTrs[part][5], log10bg);} // Estimation for relative sigma beta*gamma dependence for Ifit from TpcRS normalized to MIP
128 
129  Double_t MostProbableZShift() {return fMostProbableZShift;}
130  Double_t AverageZShift () {return fAverageZShift;}
131  Double_t I70Shift () {return fI70Shift;}
132  Double_t I60Shift () {return fI60Shift;}
133 
134  ClassDef(dEdxParameterization,0)
135 };
136 // $Id: dEdxParameterization.h,v 1.10 2015/12/24 00:16:26 fisyak Exp $
137 // $Log: dEdxParameterization.h,v $
138 // Revision 1.10 2015/12/24 00:16:26 fisyak
139 // Add TpcRS model and macros
140 //
141 #endif
dEdxParameterization(const Char_t *Tag="p10", Int_t keep3D=0, const Double_t MostProbableZShift=0, const Double_t AverageZShift=0, const Double_t I70Shift=1, const Double_t I60Shift=1)
Histograms from TpcRS dE/dx simulation for each particle type.