Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-26 02:34:37

0001 #Basic example configration file to run the b-tagging validation sequence.
0002 import FWCore.ParameterSet.Config as cms
0003 process = cms.Process("validation")
0004 
0005 # load the full reconstraction configuration, to make sure we're getting all needed dependencies
0006 process.load("Configuration.StandardSequences.MagneticField_cff")
0007 process.load("Configuration.StandardSequences.FrontierConditions_GlobalTag_cff")
0008 process.load("Configuration.StandardSequences.Reconstruction_cff")
0009 process.load("Configuration.StandardSequences.GeometryRecoDB_cff")
0010 
0011 """
0012 start customization
0013 """
0014 
0015 #Enter here the Global tags
0016 from Configuration.AlCa.GlobalTag import GlobalTag
0017 tag = GlobalTag(process.GlobalTag, 'auto:run2_mc', '')
0018 #Do you want to apply JEC? For data, no need to add 'Residual', the code is checking if events are Data or MC and add 'Residual' for Data.
0019 applyJEC = True
0020 #Data or MC?
0021 runOnMC    = True
0022 #Flavour plots for MC: "all" = plots for all jets ; "dusg" = plots for d, u, s, dus, g independently ; not mandatory and any combinations are possible 
0023 #b, c, light (dusg), non-identified (NI), PU jets plots are always produced
0024 flavPlots = "allbcldusg"
0025 #Check if jets originate from PU? option recommended (only for MC)
0026 PUid = True
0027 #List of taggers and taginfo to be considered (see example in: DQMOffline/RecoB/python/bTagCommon_cff.py)
0028 from DQMOffline.RecoB.bTagCommon_cff import *
0029 tagConfig = cms.VPSet(
0030         cms.PSet(
0031             bTagGenericAnalysisBlock,
0032             label = cms.InputTag("pfCombinedInclusiveSecondaryVertexV2BJetTags"),
0033             folder = cms.string("CSVv2")
0034         ),
0035 )
0036 
0037 """
0038 end customization
0039 """
0040 
0041 ###prints###
0042 print("is it MC ? : ", runOnMC)
0043 print("Global Tag : ", tag.globaltag)
0044 ############
0045 
0046 process.load("DQMServices.Components.DQMEnvironment_cfi")
0047 process.load("DQMServices.Core.DQM_cfg")
0048 process.load("JetMETCorrections.Configuration.JetCorrectors_cff")
0049 #keep the logging output to a nice level
0050 process.load("FWCore.MessageLogger.MessageLogger_cfi")
0051 process.MessageLogger.cerr.FwkReport.reportEvery = 100
0052 process.JECseq = cms.Sequence(process.ak4PFCHSL1FastL2L3CorrectorChain)
0053 
0054 if runOnMC:
0055     #for MC jet flavour
0056     process.load("PhysicsTools.JetMCAlgos.HadronAndPartonSelector_cfi")
0057     process.load("PhysicsTools.JetMCAlgos.AK4PFJetsMCFlavourInfos_cfi")
0058     process.ak4JetFlavourInfos.jets = cms.InputTag("ak4PFJetsCHS")
0059     process.ak4JetFlavourInfos.hadronFlavourHasPriority = cms.bool(True)
0060     process.flavourSeq = cms.Sequence(
0061         process.selectedHadronsAndPartons *
0062         process.ak4JetFlavourInfos
0063     )
0064     #Validation sequence
0065     process.load("Validation.RecoB.bTagAnalysis_cfi")
0066     process.bTagValidation.jetMCSrc = 'ak4JetFlavourInfos'
0067     process.bTagValidation.tagConfig = tagConfig
0068     process.bTagHarvestMC.tagConfig = tagConfig
0069     process.bTagValidation.flavPlots = flavPlots
0070     process.bTagHarvestMC.flavPlots = flavPlots
0071     process.bTagValidation.doPUid = cms.bool(PUid)
0072     process.bTagValidation.doJEC = applyJEC
0073     process.bTagValidation.JECsourceMC = cms.InputTag("ak4PFCHSL1FastL2L3Corrector")
0074     process.ak4GenJetsForPUid = cms.EDFilter("GenJetSelector",
0075                                              src = cms.InputTag("ak4GenJets"),
0076                                              cut = cms.string('pt > 8.'),
0077                                              filter = cms.bool(False)
0078                                              )
0079     process.load("PhysicsTools.PatAlgos.mcMatchLayer0.jetMatch_cfi")
0080     process.patJetGenJetMatch.matched = cms.InputTag("ak4GenJetsForPUid")
0081     process.patJetGenJetMatch.maxDeltaR = cms.double(0.25)
0082     process.patJetGenJetMatch.resolveAmbiguities = cms.bool(True)
0083 else :
0084     process.load("DQMOffline.RecoB.bTagAnalysisData_cfi")
0085     process.bTagAnalysis.tagConfig = tagConfig
0086     process.bTagHarvest.tagConfig = tagConfig
0087     process.bTagAnalysis.doJEC = applyJEC
0088     process.bTagAnalysis.JECsourceData = cms.InputTag("ak4PFCHSL1FastL2L3ResidualCorrector")
0089     process.JECseq *= (process.ak4PFCHSResidualCorrector * process.ak4PFCHSL1FastL2L3ResidualCorrector)
0090 
0091 
0092 process.GlobalTag = tag
0093 
0094 process.maxEvents = cms.untracked.PSet(
0095     input = cms.untracked.int32(-1)
0096 )
0097 process.source = cms.Source("PoolSource",
0098     fileNames = cms.untracked.vstring()
0099 )
0100 
0101 if runOnMC:
0102     process.dqmSeq = cms.Sequence(process.ak4GenJetsForPUid * process.patJetGenJetMatch * process.flavourSeq * process.bTagValidation * process.bTagHarvestMC * process.dqmSaver)
0103 else:
0104     process.dqmSeq = cms.Sequence(process.bTagAnalysis * process.bTagHarvest * process.dqmSaver)
0105 
0106 process.plots = cms.Path(process.JECseq*process.dqmSeq)
0107     
0108 process.dqmEnv.subSystemFolder = 'BTAG'
0109 process.dqmSaver.producer = 'DQM'
0110 process.dqmSaver.workflow = '/POG/BTAG/BJET'
0111 process.dqmSaver.convention = 'Offline'
0112 process.dqmSaver.saveByRun = cms.untracked.int32(-1)
0113 process.dqmSaver.saveAtJobEnd =cms.untracked.bool(True) 
0114 process.dqmSaver.forceRunNumber = cms.untracked.int32(1)
0115 process.PoolSource.fileNames = [
0116 
0117 ]
0118