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 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, run, ls):
58  '''
59  lumiinfo:[intlumi per LS, mean interactions ]
60 
61  intlumi is the deadtime corrected average integrated 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 
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("Run %d, LS %d: Significant probability density outside of your histogram (mean %.2f," % (run, ls, AveNumInt))
134  print("rms %.2f, integrated probability %.3f). Consider using a higher value of --maxPileupBin." % (RMSInt, totalProb))
135  else:
136  hist.Fill(AveNumInt,LSintLumi)
137  else: # have to convolute with a poisson distribution to get observed Nint
138  totalProb = 0
139  Peak = 0
140  BinWidth = hist.GetBinWidth(1)
141  for obs in range (Nbins):
142  Peak = hist.GetBinCenter(obs+1)
143  RMSWeight = ProbFromRMS[obs]
144  for bin in range (Nbins):
145  val = hist.GetBinCenter(bin+1)-0.5*BinWidth
146  prob = ROOT.TMath.Poisson (val, Peak)
147  totalProb += prob
148  hist.Fill (val, prob * LSintLumi * RMSWeight)
149 
150  if 1.0-totalProb > 0.01:
151  print("Run %d, LS %d: significant probability density outside of your histogram" % (run, ls))
152  print("Consider using a higher value of --maxPileupBin")
153 
154 
155  return hist
156 
157 
158 
159 
166 
167 if __name__ == '__main__':
168 
169  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.")
170 
171  # required
172  req_group = parser.add_argument_group('required arguments')
173  req_group.add_argument('outputfile', action='store', help='output ROOT file')
174  req_group.add_argument('-i', '--input', dest='inputfile', action='store', required=True,
175  help='input Run/LS file for your analysis in JSON format')
176  req_group.add_argument('-j', '--inputLumiJSON', dest='inputLumiJSON', action='store', required=True,
177  help='input pileup file in JSON format')
178  req_group.add_argument('-c', '--calcMode' ,dest='calcMode', action='store',
179  help='calculate either "true" or "observed" distributions',
180  choices=['true', 'observed'], required=True)
181 
182  # optional
183  parser.add_argument('-x', '--minBiasXsec', dest='minBiasXsec', action='store',
184  type=float, default=69200.0,
185  help='minimum bias cross section to use (in microbarn) (default: %(default).0f)')
186  parser.add_argument('-m', '--maxPileupBin', dest='maxPileupBin', action='store',
187  type=int, default=100, help='maximum value of pileup histogram (default: %(default)d)')
188  parser.add_argument('-n', '--numPileupBins', dest='numPileupBins', action='store',
189  type=int, default=1000, help='number of bins in pileup histogram (default: %(default)d)')
190  parser.add_argument('--pileupHistName', dest='pileupHistName', action='store',
191  default='pileup', help='name of pileup histogram (default: %(default)s)')
192  parser.add_argument('-v', '--verbose', dest='verbose', action='store_true',
193  help='verbose mode for printing' )
194 
195  options = parser.parse_args()
196  output = options.outputfile
197 
198  if options.verbose:
199  print('General configuration')
200  print('\toutputfile:' ,options.outputfile)
201  print('\tAction:' ,options.calcMode, 'luminosity distribution will be calculated')
202  print('\tinput selection file:', options.inputfile)
203  print('\tinput lumi JSON:', options.inputLumiJSON)
204  print('\tMinBiasXsec:', options.minBiasXsec)
205  print('\tmaxPileupBin:', options.maxPileupBin)
206  print('\tnumPileupBins:', options.numPileupBins)
207 
208  import ROOT
209  pileupHist = ROOT.TH1D (options.pileupHistName,options.pileupHistName,
210  options.numPileupBins,
211  0., options.maxPileupBin)
212 
213  nbins = options.numPileupBins
214  upper = options.maxPileupBin
215 
216  inpf = open(options.inputfile, 'r')
217  inputfilecontent = inpf.read()
218  inputRange = selectionParser.selectionParser (inputfilecontent).runsandls()
219 
220  #inputRange=inputFilesetParser.inputFilesetParser(options.inputfile)
221 
222  inputPileupRange=parseInputFile(options.inputLumiJSON)
223 
224  # now, we have to find the information for the input runs and lumi sections
225  # in the Lumi/Pileup list. First, loop over inputs
226 
227  for (run, lslist) in sorted (six.iteritems(inputRange)):
228  # now, look for matching run, then match lumi sections
229  # print "searching for run %d" % (run)
230  if run in inputPileupRange.keys():
231  #print run
232  LSPUlist = inputPileupRange[run]
233  # print "LSPUlist", LSPUlist
234  for LSnumber in lslist:
235  if LSnumber in LSPUlist.keys():
236  #print "found LS %d" % (LSnumber)
237  lumiInfo = LSPUlist[LSnumber]
238  # print lumiInfo
239  fillPileupHistogram(lumiInfo, options.calcMode, pileupHist,
240  options.minBiasXsec, nbins, run, LSnumber)
241  else: # trouble
242  print("Run %d, LumiSection %d not found in Lumi/Pileup input file. Check your files!" \
243  % (run,LSnumber))
244 
245  else: # trouble
246  print("Run %d not found in Lumi/Pileup input file. Check your files!" % (run))
247 
248  # print run
249  # print lslist
250 
251  histFile = ROOT.TFile.Open(output, 'recreate')
252  if not histFile:
253  raise RuntimeError("Could not open '%s' as an output root file" % output)
254  pileupHist.Write()
255  #for hist in histList:
256  # hist.Write()
257  histFile.Close()
258  print("Wrote output histogram to", output)
pileupParser.pileupParser
Definition: pileupParser.py:4
pileupCalc.MyErf
def MyErf(input)
Definition: pileupCalc.py:25
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:57
pileupCalc.parseInputFile
def parseInputFile(inputfilename)
Definition: pileupCalc.py:13
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
JetChargeProducer_cfi.exp
exp
Definition: JetChargeProducer_cfi.py:6