CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/RecoLuminosity/LumiDB/scripts/pileupCalc.py

Go to the documentation of this file.
00001 #!/usr/bin/env python
00002 VERSION='1.00'
00003 import os,sys,time
00004 import optparse
00005 from RecoLuminosity.LumiDB import pileupParser
00006 from RecoLuminosity.LumiDB import selectionParser
00007 from math import exp
00008 from math import sqrt
00009 
00010 def parseInputFile(inputfilename):
00011     '''
00012     output ({run:[ls:[inlumi, meanint]]})
00013     '''
00014     selectf=open(inputfilename,'r')
00015     inputfilecontent=selectf.read()
00016     p=pileupParser.pileupParser(inputfilecontent)                            
00017     
00018 #    p=inputFilesetParser.inputFilesetParser(inputfilename)
00019     runlsbyfile=p.runsandls()
00020     return runlsbyfile
00021 
00022 def MyErf(input):
00023 
00024     # Abramowitz and Stegun approximations for Erf (equations 7.1.25-28)
00025     X = abs(input)
00026 
00027     p = 0.47047
00028     b1 = 0.3480242
00029     b2 = -0.0958798
00030     b3 = 0.7478556
00031 
00032     T = 1.0/(1.0+p*X)
00033     cErf = 1.0 - (b1*T + b2*T*T + b3*T*T*T)*exp(-1.0*X*X)
00034     if input<0:
00035         cErf = -1.0*cErf
00036 
00037     # Alternate Erf approximation:
00038     
00039     #A1 = 0.278393
00040     #A2 = 0.230389
00041     #A3 = 0.000972
00042     #A4 = 0.078108
00043 
00044     #term = 1.0+ A1*X+ A2*X*X+ A3*X*X*X+ A4*X*X*X*X
00045     #denom = term*term*term*term
00046 
00047     #dErf = 1.0 - 1.0/denom
00048     #if input<0:
00049     #    dErf = -1.0*dErf
00050         
00051     return cErf
00052 
00053 
00054 def fillPileupHistogram (lumiInfo, calcOption, hist, minbXsec, Nbins):
00055     '''
00056     lumiinfo:[intlumi per LS, mean interactions ]
00057 
00058     intlumi is the deadtime corrected average integraged lumi per lumisection
00059     '''
00060 
00061     LSintLumi = lumiInfo[0]
00062     RMSInt = lumiInfo[1]*minbXsec
00063     AveNumInt = lumiInfo[2]*minbXsec
00064 
00065     #coeff = 0
00066 
00067     #if RMSInt > 0:
00068     #    coeff = 1.0/RMSInt/sqrt(6.283185)
00069 
00070     #expon = 2.0*RMSInt*RMSInt
00071 
00072     Sqrt2 = sqrt(2)
00073 
00074     ##Nbins = hist.GetXaxis().GetNbins()
00075 
00076     ProbFromRMS = []
00077     BinWidth = hist.GetBinWidth(1)
00078 
00079     # First, re-constitute lumi distribution for this LS from RMS:
00080     if RMSInt > 0:
00081 
00082         AreaLnew = -10.
00083         AreaL = 0
00084 
00085         for obs in range (Nbins):
00086             #Old Gaussian normalization; inaccurate for small rms and large bins
00087             #val = hist.GetBinCenter(obs+1)
00088             #prob = coeff*exp(-1.0*(val-AveNumInt)*(val-AveNumInt)/expon)
00089             #ProbFromRMS.append(prob)
00090             
00091             left = hist.GetBinLowEdge(obs+1)
00092             right = left+BinWidth
00093 
00094             argR = (AveNumInt-right)/Sqrt2/RMSInt
00095             AreaR = MyErf(argR)
00096 
00097             if AreaLnew<-5.:
00098                 argL = (AveNumInt-left)/Sqrt2/RMSInt
00099                 AreaL = MyErf(argL)
00100             else:
00101                 AreaL = AreaLnew
00102                 AreaLnew = AreaR  # save R bin value for L next time
00103 
00104             NewProb = (AreaL-AreaR)*0.5
00105 
00106             ProbFromRMS.append(NewProb)
00107 
00108             #print left, right, argL, argR, AreaL, AreaR, NewProb
00109 
00110     else:
00111         obs = hist.FindBin(AveNumInt)
00112         for bin in range (Nbins):
00113             ProbFromRMS.append(0.0)
00114         if obs<Nbins+1:            
00115             ProbFromRMS[obs] = 1.0
00116         if AveNumInt < 1.0E-5:
00117            ProbFromRMS[obs] = 0.  # just ignore zero values
00118         
00119     if calcOption == 'true':  # Just put distribution into histogram
00120         if RMSInt > 0:
00121             totalProb = 0
00122             for obs in range (Nbins):
00123                 prob = ProbFromRMS[obs]
00124                 val = hist.GetBinCenter(obs+1)
00125                 #print obs, val, RMSInt,coeff,expon,prob
00126                 totalProb += prob
00127                 hist.Fill (val, prob * LSintLumi)
00128                 
00129             if 1.0-totalProb > 0.01:
00130                 print "Significant probability density outside of your histogram"
00131                 print "Consider using a higher value of --maxPileupBin"
00132                 print "Mean %f, RMS %f, Integrated probability %f" % (AveNumInt,RMSInt,totalProb)
00133             #    hist.Fill (val, (1 - totalProb) * LSintLumi)
00134         else:
00135             hist.Fill(AveNumInt,LSintLumi)
00136     else: # have to convolute with a poisson distribution to get observed Nint
00137         totalProb = 0
00138         Peak = 0
00139         BinWidth = hist.GetBinWidth(1)
00140         for obs in range (Nbins):
00141             Peak = hist.GetBinCenter(obs+1)
00142             RMSWeight = ProbFromRMS[obs]
00143             for bin in range (Nbins):
00144                 val = hist.GetBinCenter(bin+1)-0.5*BinWidth
00145                 prob = ROOT.TMath.Poisson (val, Peak)
00146                 totalProb += prob
00147                 hist.Fill (val, prob * LSintLumi * RMSWeight)
00148 
00149         if 1.0-totalProb > 0.01:
00150             print "Significant probability density outside of your histogram"
00151             print "Consider using a higher value of --maxPileupBin"
00152 
00153 
00154     return hist
00155 
00156 
00157 
00158 ##############################
00159 ## ######################## ##
00160 ## ## ################## ## ##
00161 ## ## ## Main Program ## ## ##
00162 ## ## ################## ## ##
00163 ## ######################## ##
00164 ##############################
00165 
00166 if __name__ == '__main__':
00167 
00168     parser = optparse.OptionParser ("Usage: %prog [--options] output.root",
00169                                     description = "Script to estimate pileup distribution using xing instantaneous luminosity information and minimum bias cross section.  Output is TH1D stored in root file")
00170 #
00171 #    parser = argparse.ArgumentParser(prog=os.path.basename(sys.argv[0]),description = "Pileup Lumi Calculation",formatter_class=argparse.ArgumentDefaultsHelpFormatter)
00172     CalculationModeChoices = ['truth', 'observed']
00173 
00174     #
00175     # parse arguments
00176     #  
00177     #
00178     # basic arguments
00179     #
00180     #parser.add_argument('action',choices=allowedActions,
00181     #                    help='command actions')
00182     parser.add_option('-o',dest='outputfile',action='store',
00183                         default='PileupCalc.root',
00184                         help='output root file')
00185     parser.add_option('-i',dest='inputfile',action='store',
00186                         help='Input Run/LS file for your analysis in JSON format (required)')
00187     parser.add_option('--inputLumiJSON',dest='inputLumiJSON',action='store',
00188                         help='Input Lumi/Pileup file in JSON format (required)')
00189     parser.add_option('--calcMode',dest='calcMode',action='store',
00190                         help='Calculate either True ="true" or Observed="observed" distributions')
00191     parser.add_option('--minBiasXsec',dest='minBiasXsec',action='store',
00192                         type=float,
00193                         default=73500,
00194                         help='Minimum bias cross section assumed (in microbbarn), default %default microbarn')
00195     parser.add_option('--maxPileupBin',dest='maxPileupBin',action='store',
00196                         type=int,
00197                         default=25,
00198                         help='Maximum value of pileup histogram, default %default')
00199     parser.add_option('--numPileupBins',dest='numPileupBins',action='store',
00200                         type=int,
00201                         default=1000,
00202                         help='number of bins in pileup histogram, default %default')
00203     parser.add_option('--pileupHistName',dest='pileupHistName',action='store',
00204                         default='pileup',
00205                         help='name of pileup histogram, default %default')
00206     parser.add_option('--verbose',dest='verbose',action='store_true',help='verbose mode for printing' )
00207     
00208     # parse arguments
00209     try:
00210         (options, args) = parser.parse_args()
00211     except Exception , e:
00212         print e
00213     if not args:
00214         parser.print_usage()
00215         sys.exit()
00216     if len (args) != 1:
00217         parser.print_usage()
00218         raise RuntimeError, "Exactly one output file must be given"
00219     output = args[0]
00220     
00221 #    options=parser.parse_args()
00222 
00223     if options.verbose:
00224         print 'General configuration'
00225         print '\toutputfile: ',options.outputfile
00226         print '\tAction: ',options.calcMode, 'luminosity distribution will be calculated'
00227         print '\tinput selection file: ',options.inputfile
00228         print '\tMinBiasXsec: ',options.minBiasXsec
00229         print '\tmaxPileupBin: ',options.maxPileupBin
00230         print '\tnumPileupBins: ',options.numPileupBins
00231 
00232     import ROOT 
00233     pileupHist = ROOT.TH1D (options.pileupHistName,options.pileupHistName,
00234                       options.numPileupBins,
00235                       0., options.maxPileupBin)
00236 
00237     nbins = options.numPileupBins
00238     upper = options.maxPileupBin
00239 
00240     inpf = open (options.inputfile, 'r')
00241     inputfilecontent = inpf.read()
00242     inputRange =  selectionParser.selectionParser (inputfilecontent).runsandls()
00243 
00244     #inputRange=inputFilesetParser.inputFilesetParser(options.inputfile)
00245 
00246     if options.calcMode in ['true','observed']:
00247         inputPileupRange=parseInputFile(options.inputLumiJSON)
00248 
00249         # now, we have to find the information for the input runs and LumiSections 
00250         # in the Lumi/Pileup list. First, loop over inputs
00251 
00252         for (run, lslist) in sorted (inputRange.iteritems() ):
00253             # now, look for matching run, then match lumi sections
00254             # print "searching for run %d" % (run)
00255             if run in inputPileupRange.keys():
00256                 #print run
00257                 LSPUlist = inputPileupRange[run]
00258                 # print "LSPUlist", LSPUlist
00259                 for LSnumber in lslist:
00260                     if LSnumber in LSPUlist.keys():
00261                         #print "found LS %d" % (LSnumber)
00262                         lumiInfo = LSPUlist[LSnumber]
00263                         # print lumiInfo
00264                         fillPileupHistogram (lumiInfo, options.calcMode,
00265                                 pileupHist, options.minBiasXsec, nbins)
00266                     else: # trouble
00267                         print "Run %d, LumiSection %d not found in Lumi/Pileup input file. Check your files!" \
00268                                 % (run,LSnumber)
00269 
00270             else:  # trouble
00271                 print "Run %d not found in Lumi/Pileup input file.  Check your files!" % (run)
00272 
00273 
00274 
00275 #            print run
00276 #            print lslist
00277 
00278         histFile = ROOT.TFile.Open (output, 'recreate')
00279         if not histFile:
00280             raise RuntimeError, \
00281                  "Could not open '%s' as an output root file" % output
00282         pileupHist.Write()
00283         #for hist in histList:
00284         #    hist.Write()
00285         histFile.Close()
00286         sys.exit()
00287 
00288     else:
00289         print "must specify a pileup calculation mode via --calcMode true or --calcMode observed"
00290         sys.exit()