CMS 3D CMS Logo

applySubstructure_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 
3 from PhysicsTools.PatAlgos.tools.helpers import getPatAlgosToolsTask, addToProcessAndTask
4 
5 def applySubstructure( process, postfix="" ) :
6 
7  task = getPatAlgosToolsTask(process)
8 
9  from PhysicsTools.PatAlgos.tools.jetTools import addJetCollection
10 
11 
12  from PhysicsTools.PatAlgos.producersLayer1.jetProducer_cfi import _patJets as patJetsDefault
13 
14  #add AK8
15  addJetCollection(process, postfix=postfix, labelName = 'AK8',
16  jetSource = cms.InputTag('ak8PFJetsCHS'+postfix),
17  algo= 'AK', rParam = 0.8,
18  jetCorrections = ('AK8PFchs', cms.vstring(['L1FastJet', 'L2Relative', 'L3Absolute']), 'None'),
19  genJetCollection = cms.InputTag('slimmedGenJetsAK8')
20  )
21  getattr(process,"patJetsAK8"+postfix).userData.userFloats.src = [] # start with empty list of user floats
22  getattr(process,"selectedPatJetsAK8").cut = cms.string("pt > 170")
23 
24 
25  ## AK8 groomed masses
26  from RecoJets.Configuration.RecoPFJets_cff import ak8PFJetsCHSPruned, ak8PFJetsCHSSoftDrop
27  addToProcessAndTask('ak8PFJetsCHSPruned'+postfix, ak8PFJetsCHSPruned.clone(), process, task)
28  addToProcessAndTask('ak8PFJetsCHSSoftDrop'+postfix, ak8PFJetsCHSSoftDrop.clone(), process, task)
29  from RecoJets.JetProducers.ak8PFJetsCHS_groomingValueMaps_cfi import ak8PFJetsCHSPrunedMass, ak8PFJetsCHSTrimmedMass, ak8PFJetsCHSFilteredMass, ak8PFJetsCHSSoftDropMass
30  addToProcessAndTask('ak8PFJetsCHSPrunedMass'+postfix, ak8PFJetsCHSPrunedMass.clone(), process, task)
31  addToProcessAndTask('ak8PFJetsCHSTrimmedMass'+postfix, ak8PFJetsCHSTrimmedMass.clone(), process, task)
32  addToProcessAndTask('ak8PFJetsCHSFilteredMass'+postfix, ak8PFJetsCHSFilteredMass.clone(), process, task)
33  addToProcessAndTask('ak8PFJetsCHSSoftDropMass'+postfix, ak8PFJetsCHSSoftDropMass.clone(), process, task)
34 
35  getattr(process,"patJetsAK8").userData.userFloats.src += ['ak8PFJetsCHSPrunedMass'+postfix,'ak8PFJetsCHSSoftDropMass'+postfix]
36  getattr(process,"patJetsAK8").addTagInfos = cms.bool(False)
37 
38 
39 
40  # add Njetiness
41  process.load('RecoJets.JetProducers.nJettinessAdder_cfi')
42  task.add(process.Njettiness)
43  addToProcessAndTask('NjettinessAK8'+postfix, process.Njettiness.clone(), process, task)
44 
45 
46  getattr(process,"NjettinessAK8").src = cms.InputTag("ak8PFJetsCHS"+postfix)
47  getattr(process,"NjettinessAK8").cone = cms.double(0.8)
48  getattr(process,"patJetsAK8").userData.userFloats.src += ['NjettinessAK8'+postfix+':tau1','NjettinessAK8'+postfix+':tau2','NjettinessAK8'+postfix+':tau3']
49 
50 
51 
52 
53  #add AK8 from PUPPI
54  from RecoJets.JetProducers.ak8PFJetsPuppi_cfi import ak4PFJetsPuppi, ak8PFJetsPuppi
55  addToProcessAndTask('ak4PFJetsPuppi'+postfix,ak4PFJetsPuppi.clone(), process, task)
56  addToProcessAndTask('ak8PFJetsPuppi'+postfix,ak8PFJetsPuppi.clone(), process, task)
57  from RecoJets.Configuration.RecoPFJets_cff import ak8PFJetsPuppiSoftDrop
58  addToProcessAndTask('ak8PFJetsPuppiSoftDrop'+postfix, ak8PFJetsPuppiSoftDrop.clone(), process, task)
59  getattr(process,"ak8PFJetsPuppi").doAreaFastjet = True # even for standard ak8PFJets this is overwritten in RecoJets/Configuration/python/RecoPFJets_cff
60 
61 
62  addJetCollection(process, postfix=postfix, labelName = 'AK8Puppi',
63  jetSource = cms.InputTag('ak8PFJetsPuppi'+postfix),
64  algo= 'AK', rParam = 0.8,
65  jetCorrections = ('AK8PFPuppi', cms.vstring(['L2Relative', 'L3Absolute']), 'None'),
66  btagDiscriminators = ([x.value() for x in patJetsDefault.discriminatorSources] + ['pfBoostedDoubleSecondaryVertexAK8BJetTags']),
67  genJetCollection = cms.InputTag('slimmedGenJetsAK8')
68  )
69  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src = [] # start with empty list of user floats
70  getattr(process,"selectedPatJetsAK8Puppi"+postfix).cut = cms.string("pt > 170")
71 
72 
74  addToProcessAndTask('ak8PFJetsPuppiTracksAssociatorAtVertex'+postfix, cms.EDProducer("JetTracksAssociatorAtVertex",
75  j2tParametersVX.clone( coneSize = cms.double(0.8) ),
76  jets = cms.InputTag("ak8PFJetsPuppi") ),
77  process, task)
78  addToProcessAndTask('patJetAK8PuppiCharge'+postfix, cms.EDProducer("JetChargeProducer",
79  src = cms.InputTag("ak8PFJetsPuppiTracksAssociatorAtVertex"),
80  var = cms.string('Pt'),
81  exp = cms.double(1.0) ),
82  process, task)
83 
84  ## AK8 groomed masses
85  from RecoJets.Configuration.RecoPFJets_cff import ak8PFJetsPuppiSoftDrop
86  addToProcessAndTask('ak8PFJetsPuppiSoftDrop'+postfix, ak8PFJetsPuppiSoftDrop.clone(), process, task)
87  from RecoJets.JetProducers.ak8PFJetsPuppi_groomingValueMaps_cfi import ak8PFJetsPuppiSoftDropMass
88  addToProcessAndTask('ak8PFJetsPuppiSoftDropMass'+postfix, ak8PFJetsPuppiSoftDropMass.clone(), process, task)
89  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src += ['ak8PFJetsPuppiSoftDropMass'+postfix]
90  getattr(process,"patJetsAK8Puppi"+postfix).addTagInfos = cms.bool(False)
91 
92 
93 
94  # add Njetiness
95  addToProcessAndTask('NjettinessAK8Puppi'+postfix, process.Njettiness.clone(), process, task)
96  getattr(process,"NjettinessAK8Puppi"+postfix).src = cms.InputTag("ak8PFJetsPuppi"+postfix)
97  getattr(process,"NjettinessAK8Puppi").cone = cms.double(0.8)
98  getattr(process,"patJetsAK8Puppi").userData.userFloats.src += ['NjettinessAK8Puppi'+postfix+':tau1','NjettinessAK8Puppi'+postfix+':tau2','NjettinessAK8Puppi'+postfix+':tau3']
99 
100 
101 
102 
103  addToProcessAndTask("ak8PFJetsCHSValueMap"+postfix, cms.EDProducer("RecoJetToPatJetDeltaRValueMapProducer",
104  src = cms.InputTag("ak8PFJetsPuppi"+postfix),
105  matched = cms.InputTag("patJetsAK8"+postfix),
106  distMax = cms.double(0.8),
107  values = cms.vstring([
108  'userFloat("ak8PFJetsCHSPrunedMass"'+postfix+')',
109  'userFloat("ak8PFJetsCHSSoftDropMass"'+postfix+')',
110  'userFloat("NjettinessAK8'+postfix+':tau1")',
111  'userFloat("NjettinessAK8'+postfix+':tau2")',
112  'userFloat("NjettinessAK8'+postfix+':tau3")',
113  'pt','eta','phi','mass'
114  ]),
115  valueLabels = cms.vstring( [
116  'ak8PFJetsCHSPrunedMass',
117  'ak8PFJetsCHSSoftDropMass',
118  'NjettinessAK8CHSTau1',
119  'NjettinessAK8CHSTau2',
120  'NjettinessAK8CHSTau3',
121  'pt','eta','phi','mass'
122  ]) ),
123  process, task)
124 
125  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src += [
126  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'ak8PFJetsCHSPrunedMass'),
127  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'ak8PFJetsCHSSoftDropMass'),
128  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau1'),
129  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau2'),
130  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau3'),
131  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'pt'),
132  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'eta'),
133  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'phi'),
134  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'mass'),
135  ]
136 
137  # add Njetiness
138  process.load('RecoJets.JetProducers.nJettinessAdder_cfi')
139  task.add(process.Njettiness)
140  addToProcessAndTask('NjettinessAK8Subjets'+postfix, process.Njettiness.clone(), process, task)
141  getattr(process,"NjettinessAK8Subjets"+postfix).src = cms.InputTag("ak8PFJetsPuppiSoftDrop"+postfix, "SubJets")
142  getattr(process,"NjettinessAK8Subjets").cone = cms.double(0.8)
143 
144 
145 
146  ## PATify CHS soft drop fat jets
147  addJetCollection(
148  process,
149  postfix=postfix,
150  labelName = 'AK8PFCHSSoftDrop',
151  jetSource = cms.InputTag('ak8PFJetsCHSSoftDrop'+postfix),
152  btagDiscriminators = ['None'],
153  jetCorrections = ('AK8PFchs', ['L1FastJet', 'L2Relative', 'L3Absolute'], 'None'),
154  getJetMCFlavour = False # jet flavor disabled
155  )
156 
157 
158  ## PATify puppi soft drop fat jets
159  addJetCollection(
160  process,
161  postfix=postfix,
162  labelName = 'AK8PFPuppiSoftDrop',
163  jetSource = cms.InputTag('ak8PFJetsPuppiSoftDrop'+postfix),
164  btagDiscriminators = ['None'],
165  jetCorrections = ('AK8PFPuppi', ['L2Relative', 'L3Absolute'], 'None'),
166  getJetMCFlavour = False # jet flavor disabled
167  )
168 
169  ## PATify soft drop subjets
170  addJetCollection(
171  process,
172  postfix=postfix,
173  labelName = 'AK8PFPuppiSoftDropSubjets',
174  jetSource = cms.InputTag('ak8PFJetsPuppiSoftDrop'+postfix,'SubJets'),
175  algo = 'ak', # needed for subjet flavor clustering
176  rParam = 0.8, # needed for subjet flavor clustering
177  btagDiscriminators = ['pfDeepCSVJetTags:probb', 'pfDeepCSVJetTags:probbb', 'pfCombinedInclusiveSecondaryVertexV2BJetTags','pfCombinedMVAV2BJetTags'],
178  jetCorrections = ('AK4PFPuppi', ['L2Relative', 'L3Absolute'], 'None'),
179  explicitJTA = True, # needed for subjet b tagging
180  svClustering = True, # needed for subjet b tagging
181  genJetCollection = cms.InputTag('slimmedGenJets'),
182  fatJets=cms.InputTag('ak8PFJetsPuppi'), # needed for subjet flavor clustering
183  groomedFatJets=cms.InputTag('ak8PFJetsPuppiSoftDrop') # needed for subjet flavor clustering
184  )
185  getattr(process,"selectedPatJetsAK8PFPuppiSoftDrop"+postfix).cut = cms.string("pt > 170")
186  getattr(process,"patJetsAK8PFPuppiSoftDropSubjets"+postfix).userData.userFloats.src += ['NjettinessAK8Subjets'+postfix+':tau1','NjettinessAK8Subjets'+postfix+':tau2','NjettinessAK8Subjets'+postfix+':tau3']
187 
188  addToProcessAndTask("slimmedJetsAK8PFPuppiSoftDropSubjets"+postfix,
189  cms.EDProducer("PATJetSlimmer",
190  src = cms.InputTag("selectedPatJetsAK8PFPuppiSoftDropSubjets"),
191  packedPFCandidates = cms.InputTag("packedPFCandidates"),
192  dropJetVars = cms.string("1"),
193  dropDaughters = cms.string("0"),
194  rekeyDaughters = cms.string("1"),
195  dropTrackRefs = cms.string("1"),
196  dropSpecific = cms.string("1"),
197  dropTagInfos = cms.string("1"),
198  modifyJets = cms.bool(True),
199  mixedDaughters = cms.bool(False),
200  modifierConfig = cms.PSet( modifications = cms.VPSet() )
201  ),
202  process, task)
203 
204 
205  ## Establish references between PATified fat jets and subjets using the BoostedJetMerger
206  addToProcessAndTask("slimmedJetsAK8PFPuppiSoftDropPacked"+postfix,
207  cms.EDProducer("BoostedJetMerger",
208  jetSrc=cms.InputTag("selectedPatJetsAK8PFPuppiSoftDrop"),
209  subjetSrc=cms.InputTag("slimmedJetsAK8PFPuppiSoftDropSubjets")
210  ),
211  process, task )
212 
213 
214  addToProcessAndTask("packedPatJetsAK8"+postfix, cms.EDProducer("JetSubstructurePacker",
215  jetSrc = cms.InputTag("selectedPatJetsAK8Puppi"+postfix),
216  distMax = cms.double(0.8),
217  algoTags = cms.VInputTag(
218  cms.InputTag("slimmedJetsAK8PFPuppiSoftDropPacked"+postfix)
219  ),
220  algoLabels = cms.vstring(
221  'SoftDropPuppi'
222  ),
223  fixDaughters = cms.bool(True),
224  packedPFCandidates = cms.InputTag("packedPFCandidates"+postfix),
225  ),
226  process, task)
227 
228  # switch off daughter re-keying since it's done in the JetSubstructurePacker (and can't be done afterwards)
229  process.slimmedJetsAK8.rekeyDaughters = "0"
230 
def applySubstructure(process, postfix="")
def addToProcessAndTask(label, module, process, task)
Definition: helpers.py:27
def getPatAlgosToolsTask(process)
Definition: helpers.py:12