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 optparse
7 from RecoLuminosity.LumiDB import pileupParser
8 from RecoLuminosity.LumiDB import selectionParser
9 from math import exp
10 from math import sqrt
11 import six
12 
13 def parseInputFile(inputfilename):
14  '''
15  output ({run:[ls:[inlumi, meanint]]})
16  '''
17  selectf=open(inputfilename,'r')
18  inputfilecontent=selectf.read()
19  p=pileupParser.pileupParser(inputfilecontent)
20 
21 # p=inputFilesetParser.inputFilesetParser(inputfilename)
22  runlsbyfile=p.runsandls()
23  return runlsbyfile
24 
25 def MyErf(input):
26 
27  # Abramowitz and Stegun approximations for Erf (equations 7.1.25-28)
28  X = abs(input)
29 
30  p = 0.47047
31  b1 = 0.3480242
32  b2 = -0.0958798
33  b3 = 0.7478556
34 
35  T = 1.0/(1.0+p*X)
36  cErf = 1.0 - (b1*T + b2*T*T + b3*T*T*T)*exp(-1.0*X*X)
37  if input<0:
38  cErf = -1.0*cErf
39 
40  # Alternate Erf approximation:
41 
42  #A1 = 0.278393
43  #A2 = 0.230389
44  #A3 = 0.000972
45  #A4 = 0.078108
46 
47  #term = 1.0+ A1*X+ A2*X*X+ A3*X*X*X+ A4*X*X*X*X
48  #denom = term*term*term*term
49 
50  #dErf = 1.0 - 1.0/denom
51  #if input<0:
52  # dErf = -1.0*dErf
53 
54  return cErf
55 
56 
57 def fillPileupHistogram (lumiInfo, calcOption, hist, minbXsec, Nbins):
58  '''
59  lumiinfo:[intlumi per LS, mean interactions ]
60 
61  intlumi is the deadtime corrected average integraged lumi per lumisection
62  '''
63 
64  LSintLumi = lumiInfo[0]
65  RMSInt = lumiInfo[1]*minbXsec
66  AveNumInt = lumiInfo[2]*minbXsec
67 
68  #coeff = 0
69 
70  #if RMSInt > 0:
71  # coeff = 1.0/RMSInt/sqrt(6.283185)
72 
73  #expon = 2.0*RMSInt*RMSInt
74 
75  Sqrt2 = sqrt(2)
76 
77  ##Nbins = hist.GetXaxis().GetNbins()
78 
79  ProbFromRMS = []
80  BinWidth = hist.GetBinWidth(1)
81 
82  # First, re-constitute lumi distribution for this LS from RMS:
83  if RMSInt > 0:
84 
85  AreaLnew = -10.
86  AreaL = 0
87 
88  for obs in range (Nbins):
89  #Old Gaussian normalization; inaccurate for small rms and large bins
90  #val = hist.GetBinCenter(obs+1)
91  #prob = coeff*exp(-1.0*(val-AveNumInt)*(val-AveNumInt)/expon)
92  #ProbFromRMS.append(prob)
93 
94  left = hist.GetBinLowEdge(obs+1)
95  right = left+BinWidth
96 
97  argR = (AveNumInt-right)/Sqrt2/RMSInt
98  AreaR = MyErf(argR)
99 
100  if AreaLnew<-5.:
101  argL = (AveNumInt-left)/Sqrt2/RMSInt
102  AreaL = MyErf(argL)
103  else:
104  AreaL = AreaLnew
105  AreaLnew = AreaR # save R bin value for L next time
106 
107  NewProb = (AreaL-AreaR)*0.5
108 
109  ProbFromRMS.append(NewProb)
110 
111  #print left, right, argL, argR, AreaL, AreaR, NewProb
112 
113  else:
114  obs = hist.FindBin(AveNumInt)
115  for bin in range (Nbins):
116  ProbFromRMS.append(0.0)
117  if obs<Nbins+1:
118  ProbFromRMS[obs] = 1.0
119  if AveNumInt < 1.0E-5:
120  ProbFromRMS[obs] = 0. # just ignore zero values
121 
122  if calcOption == 'true': # Just put distribution into histogram
123  if RMSInt > 0:
124  totalProb = 0
125  for obs in range (Nbins):
126  prob = ProbFromRMS[obs]
127  val = hist.GetBinCenter(obs+1)
128  #print obs, val, RMSInt,coeff,expon,prob
129  totalProb += prob
130  hist.Fill (val, prob * LSintLumi)
131 
132  if 1.0-totalProb > 0.01:
133  print("Significant probability density outside of your histogram")
134  print("Consider using a higher value of --maxPileupBin")
135  print("Mean %f, RMS %f, Integrated probability %f" % (AveNumInt,RMSInt,totalProb))
136  # hist.Fill (val, (1 - totalProb) * LSintLumi)
137  else:
138  hist.Fill(AveNumInt,LSintLumi)
139  else: # have to convolute with a poisson distribution to get observed Nint
140  totalProb = 0
141  Peak = 0
142  BinWidth = hist.GetBinWidth(1)
143  for obs in range (Nbins):
144  Peak = hist.GetBinCenter(obs+1)
145  RMSWeight = ProbFromRMS[obs]
146  for bin in range (Nbins):
147  val = hist.GetBinCenter(bin+1)-0.5*BinWidth
148  prob = ROOT.TMath.Poisson (val, Peak)
149  totalProb += prob
150  hist.Fill (val, prob * LSintLumi * RMSWeight)
151 
152  if 1.0-totalProb > 0.01:
153  print("Significant probability density outside of your histogram")
154  print("Consider using a higher value of --maxPileupBin")
155 
156 
157  return hist
158 
159 
160 
161 ##############################
162 ## ######################## ##
163 ## ## ################## ## ##
164 ## ## ## Main Program ## ## ##
165 ## ## ################## ## ##
166 ## ######################## ##
167 ##############################
168 
169 if __name__ == '__main__':
170 
171  parser = optparse.OptionParser ("Usage: %prog [--options] output.root",
172  description = "Script to estimate pileup distribution using xing instantaneous luminosity information and minimum bias cross section. Output is TH1D stored in root file")
173 #
174 # parser = argparse.ArgumentParser(prog=os.path.basename(sys.argv[0]),description = "Pileup Lumi Calculation",formatter_class=argparse.ArgumentDefaultsHelpFormatter)
175  CalculationModeChoices = ['truth', 'observed']
176 
177  #
178  # parse arguments
179  #
180  #
181  # basic arguments
182  #
183  #parser.add_argument('action',choices=allowedActions,
184  # help='command actions')
185  parser.add_option('-o',dest='outputfile',action='store',
186  default='PileupCalc.root',
187  help='output root file')
188  parser.add_option('-i',dest='inputfile',action='store',
189  help='Input Run/LS file for your analysis in JSON format (required)')
190  parser.add_option('--inputLumiJSON',dest='inputLumiJSON',action='store',
191  help='Input Lumi/Pileup file in JSON format (required)')
192  parser.add_option('--calcMode',dest='calcMode',action='store',
193  help='Calculate either True ="true" or Observed="observed" distributions')
194  parser.add_option('--minBiasXsec',dest='minBiasXsec',action='store',
195  type=float,
196  default=73500,
197  help='Minimum bias cross section assumed (in microbbarn), default %default microbarn')
198  parser.add_option('--maxPileupBin',dest='maxPileupBin',action='store',
199  type=int,
200  default=25,
201  help='Maximum value of pileup histogram, default %default')
202  parser.add_option('--numPileupBins',dest='numPileupBins',action='store',
203  type=int,
204  default=1000,
205  help='number of bins in pileup histogram, default %default')
206  parser.add_option('--pileupHistName',dest='pileupHistName',action='store',
207  default='pileup',
208  help='name of pileup histogram, default %default')
209  parser.add_option('--verbose',dest='verbose',action='store_true',help='verbose mode for printing' )
210 
211  # parse arguments
212  try:
213  (options, args) = parser.parse_args()
214  except Exception as e:
215  print(e)
216  if not args:
217  parser.print_usage()
218  sys.exit()
219  if len (args) != 1:
220  parser.print_usage()
221  raise RuntimeError("Exactly one output file must be given")
222  output = args[0]
223 
224 # options=parser.parse_args()
225 
226  if options.verbose:
227  print('General configuration')
228  print('\toutputfile: ',options.outputfile)
229  print('\tAction: ',options.calcMode, 'luminosity distribution will be calculated')
230  print('\tinput selection file: ',options.inputfile)
231  print('\tMinBiasXsec: ',options.minBiasXsec)
232  print('\tmaxPileupBin: ',options.maxPileupBin)
233  print('\tnumPileupBins: ',options.numPileupBins)
234 
235  import ROOT
236  pileupHist = ROOT.TH1D (options.pileupHistName,options.pileupHistName,
237  options.numPileupBins,
238  0., options.maxPileupBin)
239 
240  nbins = options.numPileupBins
241  upper = options.maxPileupBin
242 
243  inpf = open (options.inputfile, 'r')
244  inputfilecontent = inpf.read()
245  inputRange = selectionParser.selectionParser (inputfilecontent).runsandls()
246 
247  #inputRange=inputFilesetParser.inputFilesetParser(options.inputfile)
248 
249  if options.calcMode in ['true','observed']:
250  inputPileupRange=parseInputFile(options.inputLumiJSON)
251 
252  # now, we have to find the information for the input runs and LumiSections
253  # in the Lumi/Pileup list. First, loop over inputs
254 
255  for (run, lslist) in sorted (six.iteritems(inputRange)):
256  # now, look for matching run, then match lumi sections
257  # print "searching for run %d" % (run)
258  if run in inputPileupRange.keys():
259  #print run
260  LSPUlist = inputPileupRange[run]
261  # print "LSPUlist", LSPUlist
262  for LSnumber in lslist:
263  if LSnumber in LSPUlist.keys():
264  #print "found LS %d" % (LSnumber)
265  lumiInfo = LSPUlist[LSnumber]
266  # print lumiInfo
267  fillPileupHistogram (lumiInfo, options.calcMode,
268  pileupHist, options.minBiasXsec, nbins)
269  else: # trouble
270  print("Run %d, LumiSection %d not found in Lumi/Pileup input file. Check your files!" \
271  % (run,LSnumber))
272 
273  else: # trouble
274  print("Run %d not found in Lumi/Pileup input file. Check your files!" % (run))
275 
276 
277 
278 # print run
279 # print lslist
280 
281  histFile = ROOT.TFile.Open (output, 'recreate')
282  if not histFile:
283  raise RuntimeError("Could not open '%s' as an output root file" % output)
284  pileupHist.Write()
285  #for hist in histList:
286  # hist.Write()
287  histFile.Close()
288  sys.exit()
289 
290  else:
291  print("must specify a pileup calculation mode via --calcMode true or --calcMode observed")
292  sys.exit()
def parseInputFile(inputfilename)
Definition: pileupCalc.py:13
def fillPileupHistogram(lumiInfo, calcOption, hist, minbXsec, Nbins)
Definition: pileupCalc.py:57
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:66
def MyErf(input)
Definition: pileupCalc.py:25
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22