CMS 3D CMS Logo

pileupCalc.py
Go to the documentation of this file.
1 #!/usr/bin/env python
2 from __future__ import print_function
3 from builtins import range
4 VERSION='1.00'
5 import os, sys, time
6 import argparse
7 from RecoLuminosity.LumiDB import pileupParser
8 from RecoLuminosity.LumiDB import selectionParser
9 from math import exp
10 from math import sqrt
11 
12 def parseInputFile(inputfilename):
13  '''
14  output ({run:[ls:[inlumi, meanint]]})
15  '''
16  selectf=open(inputfilename,'r')
17  inputfilecontent=selectf.read()
18  p=pileupParser.pileupParser(inputfilecontent)
19 
20 # p=inputFilesetParser.inputFilesetParser(inputfilename)
21  runlsbyfile=p.runsandls()
22  return runlsbyfile
23 
24 def MyErf(input):
25 
26  # Abramowitz and Stegun approximations for Erf (equations 7.1.25-28)
27  X = abs(input)
28 
29  p = 0.47047
30  b1 = 0.3480242
31  b2 = -0.0958798
32  b3 = 0.7478556
33 
34  T = 1.0/(1.0+p*X)
35  cErf = 1.0 - (b1*T + b2*T*T + b3*T*T*T)*exp(-1.0*X*X)
36  if input<0:
37  cErf = -1.0*cErf
38 
39  # Alternate Erf approximation:
40 
41  #A1 = 0.278393
42  #A2 = 0.230389
43  #A3 = 0.000972
44  #A4 = 0.078108
45 
46  #term = 1.0+ A1*X+ A2*X*X+ A3*X*X*X+ A4*X*X*X*X
47  #denom = term*term*term*term
48 
49  #dErf = 1.0 - 1.0/denom
50  #if input<0:
51  # dErf = -1.0*dErf
52 
53  return cErf
54 
55 
56 def fillPileupHistogram (lumiInfo, calcOption, hist, minbXsec, Nbins, run, ls):
57  '''
58  lumiinfo:[intlumi per LS, mean interactions ]
59 
60  intlumi is the deadtime corrected average integrated lumi per lumisection
61  '''
62 
63  LSintLumi = lumiInfo[0]
64  RMSInt = lumiInfo[1]*minbXsec
65  AveNumInt = lumiInfo[2]*minbXsec
66 
67  #coeff = 0
68 
69  #if RMSInt > 0:
70  # coeff = 1.0/RMSInt/sqrt(6.283185)
71 
72  #expon = 2.0*RMSInt*RMSInt
73 
74  Sqrt2 = sqrt(2)
75 
76 
77 
78  ProbFromRMS = []
79  BinWidth = hist.GetBinWidth(1)
80 
81  # First, re-constitute lumi distribution for this LS from RMS:
82  if RMSInt > 0:
83 
84  AreaLnew = -10.
85  AreaL = 0
86 
87  for obs in range (Nbins):
88  #Old Gaussian normalization; inaccurate for small rms and large bins
89  #val = hist.GetBinCenter(obs+1)
90  #prob = coeff*exp(-1.0*(val-AveNumInt)*(val-AveNumInt)/expon)
91  #ProbFromRMS.append(prob)
92 
93  left = hist.GetBinLowEdge(obs+1)
94  right = left+BinWidth
95 
96  argR = (AveNumInt-right)/Sqrt2/RMSInt
97  AreaR = MyErf(argR)
98 
99  if AreaLnew<-5.:
100  argL = (AveNumInt-left)/Sqrt2/RMSInt
101  AreaL = MyErf(argL)
102  else:
103  AreaL = AreaLnew
104  AreaLnew = AreaR # save R bin value for L next time
105 
106  NewProb = (AreaL-AreaR)*0.5
107 
108  ProbFromRMS.append(NewProb)
109 
110  #print left, right, argL, argR, AreaL, AreaR, NewProb
111 
112  else:
113  obs = hist.FindBin(AveNumInt)
114  for bin in range (Nbins):
115  ProbFromRMS.append(0.0)
116  if obs<Nbins+1:
117  ProbFromRMS[obs] = 1.0
118  if AveNumInt < 1.0E-5:
119  ProbFromRMS[obs] = 0. # just ignore zero values
120 
121  if calcOption == 'true': # Just put distribution into histogram
122  if RMSInt > 0:
123  totalProb = 0
124  for obs in range (Nbins):
125  prob = ProbFromRMS[obs]
126  val = hist.GetBinCenter(obs+1)
127  #print obs, val, RMSInt,coeff,expon,prob
128  totalProb += prob
129  hist.Fill (val, prob * LSintLumi)
130 
131  if 1.0-totalProb > 0.01:
132  print("Run %d, LS %d: Significant probability density outside of your histogram (mean %.2f," % (run, ls, AveNumInt))
133  print("rms %.2f, integrated probability %.3f). Consider using a higher value of --maxPileupBin." % (RMSInt, totalProb))
134  else:
135  hist.Fill(AveNumInt,LSintLumi)
136  else: # have to convolute with a poisson distribution to get observed Nint
137  totalProb = 0
138  Peak = 0
139  BinWidth = hist.GetBinWidth(1)
140  for obs in range (Nbins):
141  Peak = hist.GetBinCenter(obs+1)
142  RMSWeight = ProbFromRMS[obs]
143  for bin in range (Nbins):
144  val = hist.GetBinCenter(bin+1)-0.5*BinWidth
145  prob = ROOT.TMath.Poisson (val, Peak)
146  totalProb += prob
147  hist.Fill (val, prob * LSintLumi * RMSWeight)
148 
149  if 1.0-totalProb > 0.01:
150  print("Run %d, LS %d: significant probability density outside of your histogram" % (run, ls))
151  print("Consider using a higher value of --maxPileupBin")
152 
153 
154  return hist
155 
156 
157 
158 
165 
166 if __name__ == '__main__':
167 
168  parser = argparse.ArgumentParser(description = "Script to estimate pileup distribution using bunch instantaneous luminosity information and minimum bias cross section. Outputs a TH1D of the pileup distribution stored in a ROOT file.")
169 
170  # required
171  req_group = parser.add_argument_group('required arguments')
172  req_group.add_argument('outputfile', action='store', help='output ROOT file')
173  req_group.add_argument('-i', '--input', dest='inputfile', action='store', required=True,
174  help='input Run/LS file for your analysis in JSON format')
175  req_group.add_argument('-j', '--inputLumiJSON', dest='inputLumiJSON', action='store', required=True,
176  help='input pileup file in JSON format')
177  req_group.add_argument('-c', '--calcMode' ,dest='calcMode', action='store',
178  help='calculate either "true" or "observed" distributions',
179  choices=['true', 'observed'], required=True)
180 
181  # optional
182  parser.add_argument('-x', '--minBiasXsec', dest='minBiasXsec', action='store',
183  type=float, default=69200.0,
184  help='minimum bias cross section to use (in microbarn) (default: %(default).0f)')
185  parser.add_argument('-m', '--maxPileupBin', dest='maxPileupBin', action='store',
186  type=int, default=100, help='maximum value of pileup histogram (default: %(default)d)')
187  parser.add_argument('-n', '--numPileupBins', dest='numPileupBins', action='store',
188  type=int, default=1000, help='number of bins in pileup histogram (default: %(default)d)')
189  parser.add_argument('--pileupHistName', dest='pileupHistName', action='store',
190  default='pileup', help='name of pileup histogram (default: %(default)s)')
191  parser.add_argument('-v', '--verbose', dest='verbose', action='store_true',
192  help='verbose mode for printing' )
193 
194  options = parser.parse_args()
195  output = options.outputfile
196 
197  if options.verbose:
198  print('General configuration')
199  print('\toutputfile:' ,options.outputfile)
200  print('\tAction:' ,options.calcMode, 'luminosity distribution will be calculated')
201  print('\tinput selection file:', options.inputfile)
202  print('\tinput lumi JSON:', options.inputLumiJSON)
203  print('\tMinBiasXsec:', options.minBiasXsec)
204  print('\tmaxPileupBin:', options.maxPileupBin)
205  print('\tnumPileupBins:', options.numPileupBins)
206 
207  import ROOT
208  pileupHist = ROOT.TH1D (options.pileupHistName,options.pileupHistName,
209  options.numPileupBins,
210  0., options.maxPileupBin)
211 
212  nbins = options.numPileupBins
213  upper = options.maxPileupBin
214 
215  inpf = open(options.inputfile, 'r')
216  inputfilecontent = inpf.read()
217  inputRange = selectionParser.selectionParser (inputfilecontent).runsandls()
218 
219  #inputRange=inputFilesetParser.inputFilesetParser(options.inputfile)
220 
221  inputPileupRange=parseInputFile(options.inputLumiJSON)
222 
223  # now, we have to find the information for the input runs and lumi sections
224  # in the Lumi/Pileup list. First, loop over inputs
225 
226  for (run, lslist) in sorted (inputRange.items()):
227  # now, look for matching run, then match lumi sections
228  # print "searching for run %d" % (run)
229  if run in inputPileupRange.keys():
230  #print run
231  LSPUlist = inputPileupRange[run]
232  # print "LSPUlist", LSPUlist
233  for LSnumber in lslist:
234  if LSnumber in LSPUlist.keys():
235  #print "found LS %d" % (LSnumber)
236  lumiInfo = LSPUlist[LSnumber]
237  # print lumiInfo
238  fillPileupHistogram(lumiInfo, options.calcMode, pileupHist,
239  options.minBiasXsec, nbins, run, LSnumber)
240  else: # trouble
241  print("Run %d, LumiSection %d not found in Lumi/Pileup input file. Check your files!" \
242  % (run,LSnumber))
243 
244  else: # trouble
245  print("Run %d not found in Lumi/Pileup input file. Check your files!" % (run))
246 
247  # print run
248  # print lslist
249 
250  histFile = ROOT.TFile.Open(output, 'recreate')
251  if not histFile:
252  raise RuntimeError("Could not open '%s' as an output root file" % output)
253  pileupHist.Write()
254  #for hist in histList:
255  # hist.Write()
256  histFile.Close()
257  print("Wrote output histogram to", output)
pileupParser.pileupParser
Definition: pileupParser.py:4
pileupCalc.MyErf
def MyErf(input)
Definition: pileupCalc.py:24
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
selectionParser.selectionParser
Definition: selectionParser.py:4
print
void print(TMatrixD &m, const char *label=nullptr, bool mathematicaFormat=false)
Definition: Utilities.cc:46
pileupCalc.fillPileupHistogram
def fillPileupHistogram(lumiInfo, calcOption, hist, minbXsec, Nbins, run, ls)
Definition: pileupCalc.py:56
pileupCalc.parseInputFile
def parseInputFile(inputfilename)
Definition: pileupCalc.py:12
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6