StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StRefMultCorr.h
1 
30 #ifndef __StRefMultCorr_h__
31 #define __StRefMultCorr_h__
32 
33 // C++ headers
34 #include <vector>
35 #include <map>
36 
37 // ROOT headers
38 #include "Rtypes.h"
39 #include "TString.h"
40 #include "BadRun.h"
41 
42 //_________________
43 // Class to correct z-vertex dependence, luminosity dependence of multiplicity
45  public:
62  StRefMultCorr(const TString name="refmult", const TString subname="Def", const TString libname="Def");
64  virtual ~StRefMultCorr();
65 
67  Bool_t isBadRun(const Int_t RunId) ;
68 
69  // Event-by-event initialization. Call this function event-by-event
70  // * Default ZDC coincidence rate = 0 to make the function backward compatible
71  // --> i.e. no correction will be applied unless users set the values for 200 GeV
72  void initEvent(const UShort_t RefMult, const Double_t z, const Double_t zdcCoincidenceRate=0.0) ; // Set multiplicity, vz and zdc coincidence rate
73 
75  Bool_t isPileUpEvent(Double_t refmult, Double_t ntofmatch, Double_t vz=0., Double_t totnMIP=-999.) const;
77  Bool_t passnTofMatchRefmultCut(Double_t refmult, Double_t ntofmatch, Double_t vz=0.) const;
79  Bool_t passnTofMatchTotnMIPCut(Double_t totnMIP, Double_t ntofmatch, Double_t vz=0.) const;
80 
82  Double_t getRefMultCorr() const;
83  // Corrected multiplity, flag=0: Luminosity only, flag=1: z-vertex only, flag=2: full correction (default)
84  Double_t getRefMultCorr(const UShort_t RefMult, const Double_t z,
85  const Double_t zdcCoincidenceRate, const UInt_t flag = 2) const;
86 
88  Double_t luminosityCorrection(Double_t zdcCoincidenceRate) const;
90  Double_t vzCorrection(Double_t z) const;
92  Double_t sampleRefMult(Int_t refMult) const;
94  Double_t getShapeWeight_SubVz2Center() const;
96  Double_t triggerWeight() const;
97 
99  Double_t getWeight() const;
100 
102  Int_t getCentralityBin16() const;
104  Int_t getCentralityBin9() const;
105 
107  void init(const Int_t RunId);
108 
109  // Read scale factor from text file
110  void setVzForWeight(const Int_t nbin, const Double_t min, const Double_t max) ;
111  void readScaleForWeight(const Char_t* input) ;
112 
114  Int_t getBeginRun(const Double_t energy, const Int_t year);
116  Int_t getEndRun( const Double_t energy, const Int_t year);
117 
119  void print(const Option_t* option="") const ;
120 
122  void setVerbose(const Bool_t& verbose) { mVerbose = verbose; }
123 
124  private:
126  const TString mName;
128  const TString mSubName;
130  const TString mLibName;
131 
132  // Functions
133 
135  void clear();
137  Bool_t isIndexOk() const;
139  Bool_t isZvertexOk() const;
141  Bool_t isRefMultOk() const;
143  Bool_t isCentralityOk(const Int_t icent) const;
145  Int_t setParameterIndex(const Int_t RunId);
147  Bool_t mIsZr;
148  Bool_t mIsRu;
149 
150  // Special scale factor for Run14 to take into account the weight
151  // between different triggers
152  // - return 1 for all the other runs
153  Double_t getScaleForWeight() const ;
154 
155  // Data members
156  enum
157  {
159  mNCentrality = 16,
161  mNPar_z_vertex = 8,
163  mNPar_weight = 8,
165  mNPar_luminosity = 2
166  };
167 
168  // Use these variables to avoid varying the corrected multiplicity
169  // in the same event by random numbers
170  UShort_t mRefMult ;
171  Double_t mVz ;
172  Double_t mZdcCoincidenceRate ;
173  Double_t mRefMult_corr;
174 
175  std::vector<Int_t> mYear ;
176  std::vector<Int_t> mStart_runId ;
177  std::vector<Int_t> mStop_runId ;
178  std::vector<Double_t> mStart_zvertex ;
179  std::vector<Double_t> mStop_zvertex ;
180  std::vector<Double_t> mNormalize_stop ;
181  std::vector<Int_t> mCentrality_bins[mNCentrality+1] ;
182  std::vector<Double_t> mPar_z_vertex[mNPar_z_vertex] ;
183  std::vector<Double_t> mPar_weight[mNPar_weight] ;
184  std::vector<Double_t> mPar_luminosity[mNPar_luminosity] ;
185  Int_t mParameterIndex;
186 
187  std::multimap<std::pair<Double_t, Int_t>, Int_t> mBeginRun ;
188  std::multimap<std::pair<Double_t, Int_t>, Int_t> mEndRun ;
189  std::vector<Int_t> mBadRun ;
190 
191  // [6][680];
192  Int_t mnVzBinForWeight ;
193  std::vector<Double_t> mVzEdgeForWeight ;
194  std::vector<Double_t> mgRefMultTriggerCorrDiffVzScaleRatio ;
195 
205  Short_t mRefX;
207  Bool_t mVerbose;
208 
211  const Int_t getRefX() const;
212 
213  const Int_t getNumberOfDatasets() const;
214  void readHeaderFile();
215  void readBadRunsFromHeaderFile();
216  Int_t getVzWindowForVzDepCentDef() const;
218  Int_t getCentralityBin9VzDep() const;
219  Int_t getCentralityBin16VzDep() const;
220  std::vector<std::string> StringSplit( const std::string str, const char sep ) const;
221 
222  // Read scale factor from header file
223  void readScaleForWeight(const Int_t nRefmultBin, const Double_t *weight);
224 
227  Double_t calcPileUpRefMult(Double_t ntofmatch, Double_t x0, Double_t x1,
228  Double_t x2, Double_t x3, Double_t x4) const;
230  Bool_t isInPileUpRefMultLimits(Double_t refMult, Double_t low, Double_t hi) const
231  { return ( low < refMult && refMult < hi); }
232 
233  ClassDef(StRefMultCorr, 0)
234 };
235 #endif
236 
Double_t luminosityCorrection(Double_t zdcCoincidenceRate) const
Luminosity correction factor.
Int_t getCentralityBin9() const
Get 9 centrality bins (10% increment except for 0-5 and 5-10)
Double_t getShapeWeight_SubVz2Center() const
Shape reweighting of refmult: ratio of refMult in each Vz bin to that in the center (|Vz|&lt;10cm) ...
Int_t getCentralityBin16() const
Get 16 centrality bins (5% increment, 0-5, 5-10, ..., 75-80)
Int_t getBeginRun(const Double_t energy, const Int_t year)
Return the first runId from energy and year.
Double_t sampleRefMult(Int_t refMult) const
Sample refMult -&gt; convert integer to double.
Double_t getRefMultCorr() const
Get corrected multiplicity, correction as a function of primary z-vertex.
void print(const Option_t *option="") const
Print all parameters.
Double_t getWeight() const
Total weighting factor: incorporates shape and trigger efficiency weights.
Double_t triggerWeight() const
Trigger efficiency: fit of the Glauber/Data.
Bool_t passnTofMatchRefmultCut(Double_t refmult, Double_t ntofmatch, Double_t vz=0.) const
Check if NOT pile-up event.
Double_t vzCorrection(Double_t z) const
Vz correction factor.
Bool_t isBadRun(const Int_t RunId)
Check if run is bad.
Int_t getEndRun(const Double_t energy, const Int_t year)
Return the last runId from energy and year.
StRefMultCorr(const TString name="refmult", const TString subname="Def", const TString libname="Def")
virtual ~StRefMultCorr()
Destructor.
Bool_t isPileUpEvent(Double_t refmult, Double_t ntofmatch, Double_t vz=0., Double_t totnMIP=-999.) const
Check if pile-up event.
void setVerbose(const Bool_t &verbose)
Print debug information.
void init(const Int_t RunId)
Initialization of centrality bins etc.
Bool_t passnTofMatchTotnMIPCut(Double_t totnMIP, Double_t ntofmatch, Double_t vz=0.) const
Check if NOT pile-up event using EPD&#39;s total nMIP versus nBTOFMatch.