Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:59

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 # Define the CMSSW process
0004 process = cms.Process("RERUN")
0005 
0006 import PhysicsTools.PatAlgos.tools.helpers as configtools
0007 patAlgosToolsTask = configtools.getPatAlgosToolsTask(process)
0008 
0009 # Load the standard set of configuration modules
0010 process.load('Configuration.StandardSequences.Services_cff')
0011 patAlgosToolsTask.add(process.randomEngineStateProducer)
0012 process.load('Configuration.StandardSequences.GeometryDB_cff')
0013 process.load('Configuration.StandardSequences.MagneticField_38T_cff')
0014 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0015 
0016 # Message Logger settings
0017 process.load("FWCore.MessageService.MessageLogger_cfi")
0018 process.MessageLogger.cerr.FwkReport.reportEvery = 1
0019 
0020 # Set the process options -- Display summary at the end, enable unscheduled execution
0021 process.options = cms.untracked.PSet( 
0022     wantSummary = cms.untracked.bool(False) 
0023 )
0024 
0025 # How many events to process
0026 process.maxEvents = cms.untracked.PSet( 
0027    input = cms.untracked.int32(100)
0028 )
0029 
0030 #configurable options =======================================================================
0031 runOnData=False #data/MC switch
0032 usePrivateSQlite=False #use external JECs (sqlite file)
0033 useHFCandidates=True #create an additionnal NoHF slimmed MET collection if the option is set to false
0034 redoPuppi=True # rebuild puppiMET
0035 #===================================================================
0036 
0037 
0038 ### External JECs =====================================================================================================
0039 
0040 from Configuration.AlCa.autoCond import autoCond
0041 if runOnData:
0042   process.GlobalTag.globaltag = autoCond['run2_data']
0043 else:
0044   process.GlobalTag.globaltag = autoCond['run2_mc']
0045 
0046 
0047 #Summer16_25nsV1_MC.db
0048 
0049 if usePrivateSQlite:
0050     from CondCore.DBCommon.CondDBSetup_cfi import *
0051     import os
0052     if runOnData:
0053       era="Summer15_25nsV6_DATA"
0054     else:
0055       era="Summer15_25nsV6_MC"
0056 
0057     process.jec = cms.ESSource("PoolDBESSource",CondDBSetup,
0058                                connect = cms.string( "frontier://FrontierPrep/CMS_COND_PHYSICSTOOLS"),
0059                                #connect = cms.string('sqlite:'+era+'.db'),
0060                                toGet =  cms.VPSet(
0061             cms.PSet(
0062                 record = cms.string("JetCorrectionsRecord"),
0063                 tag = cms.string("JetCorrectorParametersCollection_"+era+"_AK4PF"),
0064                 label= cms.untracked.string("AK4PF")
0065                 ),
0066             cms.PSet(
0067                 record = cms.string("JetCorrectionsRecord"),
0068                 tag = cms.string("JetCorrectorParametersCollection_"+era+"_AK4PFchs"),
0069                 label= cms.untracked.string("AK4PFchs")
0070                 ),
0071             )
0072                                )
0073     process.es_prefer_jec = cms.ESPrefer("PoolDBESSource",'jec')
0074 
0075 
0076 ### =====================================================================================================
0077 # Define the input source
0078 if runOnData:
0079   fname = '/store/relval/CMSSW_8_0_20/MET/MINIAOD/80X_dataRun2_relval_Candidate_2016_09_02_10_27_40_RelVal_met2016B-v1/00000/2E6B9138-1C7A-E611-AE72-0025905A60DE.root' 
0080 else:
0081   fname = '/store/relval/CMSSW_8_0_20/RelValZMM_13/MINIAODSIM/80X_mcRun2_asymptotic_2016_TrancheIV_v4_Tr4GT_v4-v1/00000/64F9C946-C57A-E611-AA05-0CC47A74527A.root'
0082 
0083 # Define the input source
0084 process.source = cms.Source("PoolSource", 
0085     fileNames = cms.untracked.vstring([ fname ])
0086 )
0087 
0088 
0089 ### ---------------------------------------------------------------------------
0090 ### Removing the HF from the MET computation
0091 ### ---------------------------------------------------------------------------
0092 if not useHFCandidates:
0093     process.noHFCands = cms.EDFilter("CandPtrSelector",
0094                                      src=cms.InputTag("packedPFCandidates"),
0095                                      cut=cms.string("abs(pdgId)!=1 && abs(pdgId)!=2 && abs(eta)<3.0")
0096                                      )
0097     patAlgosToolsTask.add(process.noHFCands)
0098 
0099 #jets are rebuilt from those candidates by the tools, no need to do anything else
0100 ### =================================================================================
0101 
0102 from PhysicsTools.PatUtils.tools.runMETCorrectionsAndUncertainties import runMetCorAndUncFromMiniAOD
0103 
0104 #default configuration for miniAOD reprocessing, change the isData flag to run on data
0105 #for a full met computation, remove the pfCandColl input
0106 runMetCorAndUncFromMiniAOD(process,
0107                            isData=runOnData,
0108                            )
0109 
0110 if not useHFCandidates:
0111     runMetCorAndUncFromMiniAOD(process,
0112                                isData=runOnData,
0113                                pfCandColl=cms.InputTag("noHFCands"),
0114                                reclusterJets=True, #needed for NoHF
0115                                recoMetFromPFCs=True, #needed for NoHF
0116                                postfix="NoHF"
0117                                )
0118 
0119 if redoPuppi:
0120   from PhysicsTools.PatAlgos.slimming.puppiForMET_cff import makePuppiesFromMiniAOD
0121   makePuppiesFromMiniAOD( process );
0122 
0123   runMetCorAndUncFromMiniAOD(process,
0124                              isData=runOnData,
0125                              metType="Puppi",
0126                              recoMetFromPFCs=True,
0127                              reclusterJets=True,
0128                              jetFlavor="AK4PFPuppi",
0129                              postfix="Puppi"
0130                              )
0131 
0132 
0133 process.MINIAODSIMoutput = cms.OutputModule("PoolOutputModule",
0134     compressionLevel = cms.untracked.int32(4),
0135     compressionAlgorithm = cms.untracked.string('LZMA'),
0136     eventAutoFlushCompressedSize = cms.untracked.int32(15728640),
0137     outputCommands = cms.untracked.vstring( "keep *_slimmedMETs_*_*",
0138                                             "keep *_patPFMet_*_*",
0139                                             "keep *_patPFMetT1_*_*",
0140                                             "keep *_patPFMetT1JetResDown_*_*",
0141                                             "keep *_patPFMetT1JetResUp_*_*",
0142                                             "keep *_patPFMetT1Smear_*_*",
0143                                             "keep *_patPFMetT1SmearJetResDown_*_*",
0144                                             "keep *_patPFMetT1SmearJetResUp_*_*",
0145                                             "keep *_patPFMetT1Puppi_*_*",
0146                                             "keep *_slimmedMETsPuppi_*_*",
0147                                             ),
0148     fileName = cms.untracked.string('corMETMiniAOD.root'),
0149     dataset = cms.untracked.PSet(
0150         filterName = cms.untracked.string(''),
0151         dataTier = cms.untracked.string('')
0152     ),
0153     dropMetaData = cms.untracked.string('ALL'),
0154     fastCloning = cms.untracked.bool(False),
0155     overrideInputFileSplitLevels = cms.untracked.bool(True)
0156 )
0157 
0158 
0159 process.MINIAODSIMoutput_step = cms.EndPath(process.MINIAODSIMoutput, patAlgosToolsTask)