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 VERSION='1.00'
4 import os,sys,time
5 import optparse
6 from RecoLuminosity.LumiDB import pileupParser
7 from RecoLuminosity.LumiDB import selectionParser
8 from math import exp
9 from math import sqrt
10 import six
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):
57  '''
58  lumiinfo:[intlumi per LS, mean interactions ]
59 
60  intlumi is the deadtime corrected average integraged 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  ##Nbins = hist.GetXaxis().GetNbins()
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("Significant probability density outside of your histogram")
133  print("Consider using a higher value of --maxPileupBin")
134  print("Mean %f, RMS %f, Integrated probability %f" % (AveNumInt,RMSInt,totalProb))
135  # hist.Fill (val, (1 - totalProb) * LSintLumi)
136  else:
137  hist.Fill(AveNumInt,LSintLumi)
138  else: # have to convolute with a poisson distribution to get observed Nint
139  totalProb = 0
140  Peak = 0
141  BinWidth = hist.GetBinWidth(1)
142  for obs in range (Nbins):
143  Peak = hist.GetBinCenter(obs+1)
144  RMSWeight = ProbFromRMS[obs]
145  for bin in range (Nbins):
146  val = hist.GetBinCenter(bin+1)-0.5*BinWidth
147  prob = ROOT.TMath.Poisson (val, Peak)
148  totalProb += prob
149  hist.Fill (val, prob * LSintLumi * RMSWeight)
150 
151  if 1.0-totalProb > 0.01:
152  print("Significant probability density outside of your histogram")
153  print("Consider using a higher value of --maxPileupBin")
154 
155 
156  return hist
157 
158 
159 
160 ##############################
161 ## ######################## ##
162 ## ## ################## ## ##
163 ## ## ## Main Program ## ## ##
164 ## ## ################## ## ##
165 ## ######################## ##
166 ##############################
167 
168 if __name__ == '__main__':
169 
170  parser = optparse.OptionParser ("Usage: %prog [--options] output.root",
171  description = "Script to estimate pileup distribution using xing instantaneous luminosity information and minimum bias cross section. Output is TH1D stored in root file")
172 #
173 # parser = argparse.ArgumentParser(prog=os.path.basename(sys.argv[0]),description = "Pileup Lumi Calculation",formatter_class=argparse.ArgumentDefaultsHelpFormatter)
174  CalculationModeChoices = ['truth', 'observed']
175 
176  #
177  # parse arguments
178  #
179  #
180  # basic arguments
181  #
182  #parser.add_argument('action',choices=allowedActions,
183  # help='command actions')
184  parser.add_option('-o',dest='outputfile',action='store',
185  default='PileupCalc.root',
186  help='output root file')
187  parser.add_option('-i',dest='inputfile',action='store',
188  help='Input Run/LS file for your analysis in JSON format (required)')
189  parser.add_option('--inputLumiJSON',dest='inputLumiJSON',action='store',
190  help='Input Lumi/Pileup file in JSON format (required)')
191  parser.add_option('--calcMode',dest='calcMode',action='store',
192  help='Calculate either True ="true" or Observed="observed" distributions')
193  parser.add_option('--minBiasXsec',dest='minBiasXsec',action='store',
194  type=float,
195  default=73500,
196  help='Minimum bias cross section assumed (in microbbarn), default %default microbarn')
197  parser.add_option('--maxPileupBin',dest='maxPileupBin',action='store',
198  type=int,
199  default=25,
200  help='Maximum value of pileup histogram, default %default')
201  parser.add_option('--numPileupBins',dest='numPileupBins',action='store',
202  type=int,
203  default=1000,
204  help='number of bins in pileup histogram, default %default')
205  parser.add_option('--pileupHistName',dest='pileupHistName',action='store',
206  default='pileup',
207  help='name of pileup histogram, default %default')
208  parser.add_option('--verbose',dest='verbose',action='store_true',help='verbose mode for printing' )
209 
210  # parse arguments
211  try:
212  (options, args) = parser.parse_args()
213  except Exception as e:
214  print(e)
215  if not args:
216  parser.print_usage()
217  sys.exit()
218  if len (args) != 1:
219  parser.print_usage()
220  raise RuntimeError("Exactly one output file must be given")
221  output = args[0]
222 
223 # options=parser.parse_args()
224 
225  if options.verbose:
226  print('General configuration')
227  print('\toutputfile: ',options.outputfile)
228  print('\tAction: ',options.calcMode, 'luminosity distribution will be calculated')
229  print('\tinput selection file: ',options.inputfile)
230  print('\tMinBiasXsec: ',options.minBiasXsec)
231  print('\tmaxPileupBin: ',options.maxPileupBin)
232  print('\tnumPileupBins: ',options.numPileupBins)
233 
234  import ROOT
235  pileupHist = ROOT.TH1D (options.pileupHistName,options.pileupHistName,
236  options.numPileupBins,
237  0., options.maxPileupBin)
238 
239  nbins = options.numPileupBins
240  upper = options.maxPileupBin
241 
242  inpf = open (options.inputfile, 'r')
243  inputfilecontent = inpf.read()
244  inputRange = selectionParser.selectionParser (inputfilecontent).runsandls()
245 
246  #inputRange=inputFilesetParser.inputFilesetParser(options.inputfile)
247 
248  if options.calcMode in ['true','observed']:
249  inputPileupRange=parseInputFile(options.inputLumiJSON)
250 
251  # now, we have to find the information for the input runs and LumiSections
252  # in the Lumi/Pileup list. First, loop over inputs
253 
254  for (run, lslist) in sorted (six.iteritems(inputRange)):
255  # now, look for matching run, then match lumi sections
256  # print "searching for run %d" % (run)
257  if run in inputPileupRange.keys():
258  #print run
259  LSPUlist = inputPileupRange[run]
260  # print "LSPUlist", LSPUlist
261  for LSnumber in lslist:
262  if LSnumber in LSPUlist.keys():
263  #print "found LS %d" % (LSnumber)
264  lumiInfo = LSPUlist[LSnumber]
265  # print lumiInfo
266  fillPileupHistogram (lumiInfo, options.calcMode,
267  pileupHist, options.minBiasXsec, nbins)
268  else: # trouble
269  print("Run %d, LumiSection %d not found in Lumi/Pileup input file. Check your files!" \
270  % (run,LSnumber))
271 
272  else: # trouble
273  print("Run %d not found in Lumi/Pileup input file. Check your files!" % (run))
274 
275 
276 
277 # print run
278 # print lslist
279 
280  histFile = ROOT.TFile.Open (output, 'recreate')
281  if not histFile:
282  raise RuntimeError("Could not open '%s' as an output root file" % output)
283  pileupHist.Write()
284  #for hist in histList:
285  # hist.Write()
286  histFile.Close()
287  sys.exit()
288 
289  else:
290  print("must specify a pileup calculation mode via --calcMode true or --calcMode observed")
291  sys.exit()
def parseInputFile(inputfilename)
Definition: pileupCalc.py:12
def fillPileupHistogram(lumiInfo, calcOption, hist, minbXsec, Nbins)
Definition: pileupCalc.py:56
S & print(S &os, JobReport::InputFile const &f)
Definition: JobReport.cc:65
def MyErf(input)
Definition: pileupCalc.py:24
T sqrt(T t)
Definition: SSEVec.h:18
Abs< T >::type abs(const T &t)
Definition: Abs.h:22