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 
15  # Configure the RECO jets
16  from RecoJets.JetProducers.ak4PFJets_cfi import ak4PFJetsPuppi
17  from RecoJets.JetProducers.ak8PFJets_cfi import ak8PFJetsPuppi, ak8PFJetsPuppiSoftDrop, ak8PFJetsPuppiConstituents, ak8PFJetsCHSConstituents
18  from RecoJets.JetProducers.ak8GenJets_cfi import ak8GenJets, ak8GenJetsSoftDrop, ak8GenJetsConstituents
19  addToProcessAndTask('ak4PFJetsPuppi'+postfix,ak4PFJetsPuppi.clone(), process, task)
20  addToProcessAndTask('ak8PFJetsPuppi'+postfix,ak8PFJetsPuppi.clone(), process, task)
21  addToProcessAndTask('ak8PFJetsPuppiConstituents', ak8PFJetsPuppiConstituents.clone(cut = cms.string('pt > 170.0 && abs(rapidity()) < 2.4') ), process, task )
22  addToProcessAndTask('ak8PFJetsCHSConstituents', ak8PFJetsCHSConstituents.clone(), process, task )
23  addToProcessAndTask('ak8PFJetsPuppiSoftDrop'+postfix, ak8PFJetsPuppiSoftDrop.clone( src = cms.InputTag('ak8PFJetsPuppiConstituents', 'constituents') ), process, task)
24  addToProcessAndTask('ak8GenJetsNoNuConstituents'+postfix, ak8GenJetsConstituents.clone(src='ak8GenJetsNoNu'), process, task )
25  addToProcessAndTask('ak8GenJetsNoNuSoftDrop'+postfix,ak8GenJetsSoftDrop.clone(src=cms.InputTag('ak8GenJetsNoNuConstituents'+postfix, 'constituents')),process,task)
26  addToProcessAndTask('slimmedGenJetsAK8SoftDropSubJets'+postfix,
27  cms.EDProducer("PATGenJetSlimmer",
28  src = cms.InputTag("ak8GenJetsNoNuSoftDrop"+postfix, "SubJets"),
29  packedGenParticles = cms.InputTag("packedGenParticles"),
30  cut = cms.string(""),
31  cutLoose = cms.string(""),
32  nLoose = cms.uint32(0),
33  clearDaughters = cms.bool(False), #False means rekeying
34  dropSpecific = cms.bool(True), # Save space
35  ), process, task )
36 
37 
38  #add AK8 CHS
39  addJetCollection(process, postfix=postfix, labelName = 'AK8',
40  jetSource = cms.InputTag('ak8PFJetsCHS'+postfix),
41  algo= 'AK', rParam = 0.8,
42  btagDiscriminators = ['None'],
43  jetCorrections = ('AK8PFchs', cms.vstring(['L1FastJet', 'L2Relative', 'L3Absolute']), 'None'),
44  genJetCollection = cms.InputTag('slimmedGenJetsAK8')
45  )
46  getattr(process,"patJetsAK8"+postfix).userData.userFloats.src = [] # start with empty list of user floats
47  getattr(process,"selectedPatJetsAK8").cut = cms.string("pt > 170")
48 
49 
50  ## add AK8 groomed masses with CHS
51  from RecoJets.Configuration.RecoPFJets_cff import ak8PFJetsCHSPruned, ak8PFJetsCHSSoftDrop
52  addToProcessAndTask('ak8PFJetsCHSPruned'+postfix, ak8PFJetsCHSPruned.clone(), process, task)
53  addToProcessAndTask('ak8PFJetsCHSSoftDrop'+postfix, ak8PFJetsCHSSoftDrop.clone(), process, task)
54  from RecoJets.JetProducers.ak8PFJetsCHS_groomingValueMaps_cfi import ak8PFJetsCHSPrunedMass, ak8PFJetsCHSTrimmedMass, ak8PFJetsCHSFilteredMass, ak8PFJetsCHSSoftDropMass
55  addToProcessAndTask('ak8PFJetsCHSPrunedMass'+postfix, ak8PFJetsCHSPrunedMass.clone(), process, task)
56  addToProcessAndTask('ak8PFJetsCHSTrimmedMass'+postfix, ak8PFJetsCHSTrimmedMass.clone(), process, task)
57  addToProcessAndTask('ak8PFJetsCHSFilteredMass'+postfix, ak8PFJetsCHSFilteredMass.clone(), process, task)
58  addToProcessAndTask('ak8PFJetsCHSSoftDropMass'+postfix, ak8PFJetsCHSSoftDropMass.clone(), process, task)
59 
60  getattr(process,"patJetsAK8").userData.userFloats.src += ['ak8PFJetsCHSPrunedMass'+postfix,'ak8PFJetsCHSSoftDropMass'+postfix]
61  getattr(process,"patJetsAK8").addTagInfos = cms.bool(False)
62 
63  # add Njetiness for CHS
64  process.load('RecoJets.JetProducers.nJettinessAdder_cfi')
65  task.add(process.Njettiness)
66  addToProcessAndTask('NjettinessAK8'+postfix, process.Njettiness.clone(), process, task)
67  getattr(process,"NjettinessAK8").src = cms.InputTag("ak8PFJetsCHS"+postfix)
68  getattr(process,"NjettinessAK8").cone = cms.double(0.8)
69  getattr(process,"patJetsAK8").userData.userFloats.src += ['NjettinessAK8'+postfix+':tau1','NjettinessAK8'+postfix+':tau2','NjettinessAK8'+postfix+':tau3','NjettinessAK8'+postfix+':tau4']
70 
71  # add Njetiness from CHS
72  addToProcessAndTask('NjettinessAK8Subjets'+postfix, process.Njettiness.clone(), process, task)
73  getattr(process,"NjettinessAK8Subjets"+postfix).src = cms.InputTag("ak8PFJetsPuppiSoftDrop"+postfix, "SubJets")
74  getattr(process,"NjettinessAK8Subjets").cone = cms.double(0.8)
75 
76  ## PATify CHS soft drop fat jets
77  addJetCollection(
78  process,
79  postfix=postfix,
80  labelName = 'AK8PFCHSSoftDrop',
81  jetSource = cms.InputTag('ak8PFJetsCHSSoftDrop'+postfix),
82  btagDiscriminators = ['None'],
83  jetCorrections = ('AK8PFchs', ['L1FastJet', 'L2Relative', 'L3Absolute'], 'None'),
84  getJetMCFlavour = False # jet flavor disabled
85  )
86 
87 
88 
89 
90 
91  #add RECO AK8 from PUPPI and RECO AK8 PUPPI with soft drop... will be needed by ungroomed AK8 jets later
92  ## PATify puppi soft drop fat jets
93  addJetCollection(
94  process,
95  postfix=postfix,
96  labelName = 'AK8PFPuppiSoftDrop' + postfix,
97  jetSource = cms.InputTag('ak8PFJetsPuppiSoftDrop'+postfix),
98  btagDiscriminators = ['None'],
99  genJetCollection = cms.InputTag('slimmedGenJetsAK8'),
100  jetCorrections = ('AK8PFPuppi', ['L2Relative', 'L3Absolute'], 'None'),
101  getJetMCFlavour = False # jet flavor disabled
102  )
103  ## PATify soft drop subjets
104  addJetCollection(
105  process,
106  postfix=postfix,
107  labelName = 'AK8PFPuppiSoftDropSubjets',
108  jetSource = cms.InputTag('ak8PFJetsPuppiSoftDrop'+postfix,'SubJets'),
109  algo = 'ak', # needed for subjet flavor clustering
110  rParam = 0.8, # needed for subjet flavor clustering
111  btagDiscriminators = ['pfDeepCSVJetTags:probb', 'pfDeepCSVJetTags:probbb', 'pfCombinedInclusiveSecondaryVertexV2BJetTags','pfCombinedMVAV2BJetTags'],
112  jetCorrections = ('AK4PFPuppi', ['L2Relative', 'L3Absolute'], 'None'),
113  explicitJTA = True, # needed for subjet b tagging
114  svClustering = True, # needed for subjet b tagging
115  genJetCollection = cms.InputTag('slimmedGenJetsAK8SoftDropSubJets'),
116  fatJets=cms.InputTag('ak8PFJetsPuppi'), # needed for subjet flavor clustering
117  groomedFatJets=cms.InputTag('ak8PFJetsPuppiSoftDrop') # needed for subjet flavor clustering
118  )
119 
120 
121  # add groomed ECFs and N-subjettiness to soft dropped pat::Jets for fat jets and subjets
122  process.load('RecoJets.JetProducers.ECF_cff')
123  addToProcessAndTask('nb1AK8PuppiSoftDrop'+postfix, process.ecfNbeta1.clone(src = cms.InputTag("ak8PFJetsPuppiSoftDrop"+postfix), cuts = cms.vstring('', '', 'pt > 250')), process, task)
124  addToProcessAndTask('nb2AK8PuppiSoftDrop'+postfix, process.ecfNbeta2.clone(src = cms.InputTag("ak8PFJetsPuppiSoftDrop"+postfix), cuts = cms.vstring('', '', 'pt > 250')), process, task)
125 
126  #too slow now ==> disable
127  from Configuration.Eras.Modifier_pp_on_AA_2018_cff import pp_on_AA_2018
128  from Configuration.Eras.Modifier_pp_on_XeXe_2017_cff import pp_on_XeXe_2017
129  from Configuration.Eras.Modifier_phase2_common_cff import phase2_common
130  for e in [pp_on_XeXe_2017, pp_on_AA_2018, phase2_common]:
131  e.toModify(getattr(process,'nb1AK8PuppiSoftDrop'+postfix), cuts = ['pt > 999999', 'pt > 999999', 'pt > 999999'] )
132  e.toModify(getattr(process,'nb2AK8PuppiSoftDrop'+postfix), cuts = ['pt > 999999', 'pt > 999999', 'pt > 999999'] )
133 
134 
135  getattr(process,"patJetsAK8PFPuppiSoftDrop").userData.userFloats.src += ['nb1AK8PuppiSoftDrop'+postfix+':ecfN2','nb1AK8PuppiSoftDrop'+postfix+':ecfN3']
136  getattr(process,"patJetsAK8PFPuppiSoftDrop").userData.userFloats.src += ['nb2AK8PuppiSoftDrop'+postfix+':ecfN2','nb2AK8PuppiSoftDrop'+postfix+':ecfN3']
137  addToProcessAndTask('nb1AK8PuppiSoftDropSubjets'+postfix, process.ecfNbeta1.clone(src = cms.InputTag("ak8PFJetsPuppiSoftDrop"+postfix, "SubJets")), process, task)
138  addToProcessAndTask('nb2AK8PuppiSoftDropSubjets'+postfix, process.ecfNbeta2.clone(src = cms.InputTag("ak8PFJetsPuppiSoftDrop"+postfix, "SubJets")), process, task)
139  getattr(process,"patJetsAK8PFPuppiSoftDropSubjets"+postfix).userData.userFloats.src += ['nb1AK8PuppiSoftDropSubjets'+postfix+':ecfN2','nb1AK8PuppiSoftDropSubjets'+postfix+':ecfN3']
140  getattr(process,"patJetsAK8PFPuppiSoftDropSubjets"+postfix).userData.userFloats.src += ['nb2AK8PuppiSoftDropSubjets'+postfix+':ecfN2','nb2AK8PuppiSoftDropSubjets'+postfix+':ecfN3']
141  getattr(process,"patJetsAK8PFPuppiSoftDropSubjets"+postfix).userData.userFloats.src += ['NjettinessAK8Subjets'+postfix+':tau1','NjettinessAK8Subjets'+postfix+':tau2','NjettinessAK8Subjets'+postfix+':tau3','NjettinessAK8Subjets'+postfix+':tau4']
142 
143  for e in [pp_on_XeXe_2017, pp_on_AA_2018, phase2_common]:
144  e.toModify(getattr(process,'nb1AK8PuppiSoftDropSubjets'+postfix), cuts = ['pt > 999999', 'pt > 999999', 'pt > 999999'] )
145  e.toModify(getattr(process,'nb2AK8PuppiSoftDropSubjets'+postfix), cuts = ['pt > 999999', 'pt > 999999', 'pt > 999999'] )
146 
147  # rekey the groomed ECF value maps to the ungroomed reco jets, which will then be picked
148  # up by PAT in the user floats.
149  addToProcessAndTask("ak8PFJetsPuppiSoftDropValueMap"+postfix,
150  cms.EDProducer("RecoJetToPatJetDeltaRValueMapProducer",
151  src = cms.InputTag("ak8PFJetsPuppi"+postfix),
152  matched = cms.InputTag("patJetsAK8PFPuppiSoftDrop"+postfix),
153  distMax = cms.double(0.8),
154  values = cms.vstring([
155  'userFloat("nb1AK8PuppiSoftDrop'+postfix+':ecfN2")',
156  'userFloat("nb1AK8PuppiSoftDrop'+postfix+':ecfN3")',
157  'userFloat("nb2AK8PuppiSoftDrop'+postfix+':ecfN2")',
158  'userFloat("nb2AK8PuppiSoftDrop'+postfix+':ecfN3")',
159  ]),
160  valueLabels = cms.vstring( [
161  'nb1AK8PuppiSoftDropN2',
162  'nb1AK8PuppiSoftDropN3',
163  'nb2AK8PuppiSoftDropN2',
164  'nb2AK8PuppiSoftDropN3',
165  ]) ),
166  process, task)
167 
168 
169  # Patify AK8 PF PUPPI
170  addJetCollection(process, postfix=postfix, labelName = 'AK8Puppi',
171  jetSource = cms.InputTag('ak8PFJetsPuppi'+postfix),
172  algo= 'AK', rParam = 0.8,
173  jetCorrections = ('AK8PFPuppi', cms.vstring(['L2Relative', 'L3Absolute']), 'None'),
174  btagDiscriminators = ([
175  'pfCombinedSecondaryVertexV2BJetTags',
176  'pfCombinedInclusiveSecondaryVertexV2BJetTags',
177  'pfCombinedMVAV2BJetTags',
178  'pfDeepCSVJetTags:probb',
179  'pfDeepCSVJetTags:probc',
180  'pfDeepCSVJetTags:probudsg',
181  'pfDeepCSVJetTags:probbb',
182  'pfBoostedDoubleSecondaryVertexAK8BJetTags']),
183  genJetCollection = cms.InputTag('slimmedGenJetsAK8')
184  )
185  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src = [] # start with empty list of user floats
186  getattr(process,"selectedPatJetsAK8Puppi"+postfix).cut = cms.string("pt > 100")
187  getattr(process,"selectedPatJetsAK8Puppi"+postfix).cutLoose = cms.string("pt > 30")
188  getattr(process,"selectedPatJetsAK8Puppi"+postfix).nLoose = cms.uint32(3)
189 
190  from RecoJets.JetAssociationProducers.j2tParametersVX_cfi import j2tParametersVX
191  addToProcessAndTask('ak8PFJetsPuppiTracksAssociatorAtVertex'+postfix, cms.EDProducer("JetTracksAssociatorAtVertex",
192  j2tParametersVX.clone( coneSize = cms.double(0.8) ),
193  jets = cms.InputTag("ak8PFJetsPuppi") ),
194  process, task)
195  addToProcessAndTask('patJetAK8PuppiCharge'+postfix, cms.EDProducer("JetChargeProducer",
196  src = cms.InputTag("ak8PFJetsPuppiTracksAssociatorAtVertex"),
197  var = cms.string('Pt'),
198  exp = cms.double(1.0) ),
199  process, task)
200 
201  ## now add AK8 groomed masses and ECF
202  from RecoJets.JetProducers.ak8PFJetsPuppi_groomingValueMaps_cfi import ak8PFJetsPuppiSoftDropMass
203  addToProcessAndTask('ak8PFJetsPuppiSoftDropMass'+postfix, ak8PFJetsPuppiSoftDropMass.clone(), process, task)
204  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src += ['ak8PFJetsPuppiSoftDropMass'+postfix]
205  getattr(process,"patJetsAK8Puppi"+postfix).addTagInfos = cms.bool(False)
206  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src += [
207  cms.InputTag('ak8PFJetsPuppiSoftDropValueMap'+postfix,'nb1AK8PuppiSoftDropN2'),
208  cms.InputTag('ak8PFJetsPuppiSoftDropValueMap'+postfix,'nb1AK8PuppiSoftDropN3'),
209  cms.InputTag('ak8PFJetsPuppiSoftDropValueMap'+postfix,'nb2AK8PuppiSoftDropN2'),
210  cms.InputTag('ak8PFJetsPuppiSoftDropValueMap'+postfix,'nb2AK8PuppiSoftDropN3'),
211  ]
212 
213 
214  # add PUPPI Njetiness
215  addToProcessAndTask('NjettinessAK8Puppi'+postfix, process.Njettiness.clone(), process, task)
216  getattr(process,"NjettinessAK8Puppi"+postfix).src = cms.InputTag("ak8PFJetsPuppi"+postfix)
217  getattr(process,"NjettinessAK8Puppi").cone = cms.double(0.8)
218  getattr(process,"patJetsAK8Puppi").userData.userFloats.src += ['NjettinessAK8Puppi'+postfix+':tau1','NjettinessAK8Puppi'+postfix+':tau2','NjettinessAK8Puppi'+postfix+':tau3','NjettinessAK8Puppi'+postfix+':tau4']
219 
220  # Now combine the CHS and PUPPI information into the PUPPI jets via delta R value maps
221  addToProcessAndTask("ak8PFJetsCHSValueMap"+postfix, cms.EDProducer("RecoJetToPatJetDeltaRValueMapProducer",
222  src = cms.InputTag("ak8PFJetsPuppi"+postfix),
223  matched = cms.InputTag("patJetsAK8"+postfix),
224  distMax = cms.double(0.8),
225  values = cms.vstring([
226  'userFloat("ak8PFJetsCHSPrunedMass"'+postfix+')',
227  'userFloat("ak8PFJetsCHSSoftDropMass"'+postfix+')',
228  'userFloat("NjettinessAK8'+postfix+':tau1")',
229  'userFloat("NjettinessAK8'+postfix+':tau2")',
230  'userFloat("NjettinessAK8'+postfix+':tau3")',
231  'userFloat("NjettinessAK8'+postfix+':tau4")',
232  'pt','eta','phi','mass', 'jetArea', 'jecFactor(0)'
233  ]),
234  valueLabels = cms.vstring( [
235  'ak8PFJetsCHSPrunedMass',
236  'ak8PFJetsCHSSoftDropMass',
237  'NjettinessAK8CHSTau1',
238  'NjettinessAK8CHSTau2',
239  'NjettinessAK8CHSTau3',
240  'NjettinessAK8CHSTau4',
241  'pt','eta','phi','mass', 'jetArea', 'rawFactor'
242  ]) ),
243  process, task)
244 
245 
246  # Now set up the user floats
247  getattr(process,"patJetsAK8Puppi"+postfix).userData.userFloats.src += [
248  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'ak8PFJetsCHSPrunedMass'),
249  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'ak8PFJetsCHSSoftDropMass'),
250  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau1'),
251  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau2'),
252  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau3'),
253  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'NjettinessAK8CHSTau4'),
254  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'pt'),
255  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'eta'),
256  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'phi'),
257  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'mass'),
258  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'jetArea'),
259  cms.InputTag('ak8PFJetsCHSValueMap'+postfix,'rawFactor'),
260  ]
261 
262 
263 
264 
265  addToProcessAndTask("slimmedJetsAK8PFPuppiSoftDropSubjets"+postfix,
266  cms.EDProducer("PATJetSlimmer",
267  src = cms.InputTag("selectedPatJetsAK8PFPuppiSoftDropSubjets"),
268  packedPFCandidates = cms.InputTag("packedPFCandidates"),
269  dropJetVars = cms.string("1"),
270  dropDaughters = cms.string("0"),
271  rekeyDaughters = cms.string("1"),
272  dropTrackRefs = cms.string("1"),
273  dropSpecific = cms.string("1"),
274  dropTagInfos = cms.string("1"),
275  modifyJets = cms.bool(True),
276  mixedDaughters = cms.bool(False),
277  modifierConfig = cms.PSet( modifications = cms.VPSet() )
278  ),
279  process, task)
280 
281 
282  ## Establish references between PATified fat jets and subjets using the BoostedJetMerger
283  addToProcessAndTask("slimmedJetsAK8PFPuppiSoftDropPacked"+postfix,
284  cms.EDProducer("BoostedJetMerger",
285  jetSrc=cms.InputTag("selectedPatJetsAK8PFPuppiSoftDrop"),
286  subjetSrc=cms.InputTag("slimmedJetsAK8PFPuppiSoftDropSubjets")
287  ),
288  process, task )
289 
290 
291  addToProcessAndTask("packedPatJetsAK8"+postfix, cms.EDProducer("JetSubstructurePacker",
292  jetSrc = cms.InputTag("selectedPatJetsAK8Puppi"+postfix),
293  distMax = cms.double(0.8),
294  algoTags = cms.VInputTag(
295  cms.InputTag("slimmedJetsAK8PFPuppiSoftDropPacked"+postfix)
296  ),
297  algoLabels = cms.vstring(
298  'SoftDropPuppi'
299  ),
300  fixDaughters = cms.bool(True),
301  packedPFCandidates = cms.InputTag("packedPFCandidates"+postfix),
302  ),
303  process, task)
304 
305  # switch off daughter re-keying since it's done in the JetSubstructurePacker (and can't be done afterwards)
306  process.slimmedJetsAK8.rekeyDaughters = "0"
307  # Reconfigure the slimmedAK8 jet information to keep
308  process.slimmedJetsAK8.dropDaughters = cms.string("pt < 170")
309  process.slimmedJetsAK8.dropSpecific = cms.string("pt < 170")
310  process.slimmedJetsAK8.dropTagInfos = cms.string("pt < 170")
def applySubstructure(process, postfix="")
def addToProcessAndTask(label, module, process, task)
Definition: helpers.py:29
def getPatAlgosToolsTask(process)
Definition: helpers.py:14