CMS 3D CMS Logo

nanogen_cff.py
Go to the documentation of this file.
3 from PhysicsTools.NanoAOD.globals_cff import genTable
4 from PhysicsTools.NanoAOD.met_cff import metMCTable
9 from PhysicsTools.NanoAOD.common_cff import Var,CandVars
11 from PhysicsTools.NanoAOD.common_cff import CandVars
12 
13 nanoMetadata = cms.EDProducer("UniqueStringProducer",
14  strings = cms.PSet(
15  tag = cms.string("untagged"),
16  )
17 )
18 
19 nanogenSequence = cms.Sequence(
20  nanoMetadata+
21  particleLevel+
22  genJetTable+
23  patJetPartons+
24  genJetFlavourAssociation+
25  genJetFlavourTable+
26  genJetAK8Table+
27  genJetAK8FlavourAssociation+
28  genJetAK8FlavourTable+
29  tauGenJets+
30  tauGenJetsSelectorAllHadrons+
31  genVisTaus+
32  genVisTauTable+
33  genTable+
34  genParticleTables+
35  genVertexTables+
36  tautagger+
37  rivetProducerHTXS+
38  particleLevelTables+
39  metMCTable+
40  genWeightsTable+
41  lheInfoTable
42 )
43 
45  process.rivetMetTable.extension = False
46  process.lheInfoTable.storeLHEParticles = True
47  process.lheInfoTable.precision = 14
48  process.genWeightsTable.keepAllPSWeights = True
49  process.genJetFlavourAssociation.jets = process.genJetTable.src
50  process.genJetFlavourTable.src = process.genJetTable.src
51  process.genJetAK8FlavourAssociation.jets = process.genJetAK8Table.src
52  process.genJetAK8FlavourTable.src = process.genJetAK8Table.src
53  process.particleLevel.particleMaxEta = 999.
54  process.particleLevel.lepMinPt = 0.
55  process.particleLevel.lepMaxEta = 999.
56  process.genJetFlavourTable.jetFlavourInfos = "genJetFlavourAssociation"
57  # Same as default RECO
58  setGenPtPrecision(process, CandVars.pt.precision)
59  setGenEtaPrecision(process, CandVars.eta.precision)
60  setGenPhiPrecision(process, CandVars.phi.precision)
61 
63  process.nanogenSequence.insert(0, process.genParticles2HepMCHiggsVtx)
64  process.nanogenSequence.insert(0, process.genParticles2HepMC)
65  process.nanogenSequence.insert(0, process.mergedGenParticles)
66 
67  from Configuration.Eras.Modifier_run2_miniAOD_80XLegacy_cff import run2_miniAOD_80XLegacy
68  from Configuration.Eras.Modifier_run2_nanoAOD_94X2016_cff import run2_nanoAOD_94X2016
69  from Configuration.Eras.Modifier_run2_nanoAOD_94XMiniAODv1_cff import run2_nanoAOD_94XMiniAODv1
70  from Configuration.Eras.Modifier_run2_nanoAOD_94XMiniAODv2_cff import run2_nanoAOD_94XMiniAODv2
71  from Configuration.Eras.Modifier_run2_nanoAOD_102Xv1_cff import run2_nanoAOD_102Xv1
72  from Configuration.Eras.Modifier_run2_nanoAOD_92X_cff import run2_nanoAOD_92X
73 
74  (run2_nanoAOD_92X | run2_miniAOD_80XLegacy | run2_nanoAOD_94X2016 | run2_nanoAOD_94X2016 | \
75  run2_nanoAOD_94XMiniAODv1 | run2_nanoAOD_94XMiniAODv2 | \
76  run2_nanoAOD_102Xv1).toReplaceWith(nanogenSequence, nanogenSequence.copyAndExclude([genVertexTable, genVertexT0Table]))
77 
78  process.metMCTable.src = "slimmedMETs"
79  process.metMCTable.variables.pt = Var("genMET.pt", float, doc="pt")
80  process.metMCTable.variables.phi = Var("genMET.phi", float, doc="phi")
81  process.metMCTable.variables.phi.precision = CandVars.phi.precision
82 
83  process.rivetProducerHTXS.HepMCCollection = "genParticles2HepMCHiggsVtx:unsmeared"
84  process.genParticleTable.src = "prunedGenParticles"
85  process.patJetPartons.particles = "prunedGenParticles"
86  process.particleLevel.src = "genParticles2HepMC:unsmeared"
87 
88  process.genJetTable.src = "slimmedGenJets"
89  process.genJetAK8Table.src = "slimmedGenJetsAK8"
90  process.tauGenJets.GenParticles = "prunedGenParticles"
91  process.genVisTaus.srcGenParticles = "prunedGenParticles"
92 
93  nanoGenCommonCustomize(process)
94 
95  return process
96 
97 def customizeNanoGEN(process):
98  process.metMCTable.src = "genMetTrue"
99  process.metMCTable.variables = cms.PSet(PTVars)
100 
101  process.rivetProducerHTXS.HepMCCollection = "generatorSmeared"
102  process.genParticleTable.src = "genParticles"
103  process.patJetPartons.particles = "genParticles"
104  process.particleLevel.src = "generatorSmeared"
105 
106  process.genJetTable.src = "ak4GenJets"
107  process.genJetAK8Table.src = "ak8GenJets"
108  process.tauGenJets.GenParticles = "genParticles"
109  process.genVisTaus.srcGenParticles = "genParticles"
110 
111  # In case customizeNanoGENFromMini has already been called
112  process.nanogenSequence.remove(process.genParticles2HepMCHiggsVtx)
113  process.nanogenSequence.remove(process.genParticles2HepMC)
114  process.nanogenSequence.remove(process.mergedGenParticles)
115  nanoGenCommonCustomize(process)
116  return process
117 
118 # Prune gen particles with tight conditions applied in usual NanoAOD
120  process.finalGenParticles = finalGenParticles.clone()
121  process.genParticleTable.src = "prunedGenParticles"
122  process.patJetPartons.particles = "prunedGenParticles"
123  process.nanogenSequence.insert(0, process.finalGenParticles)
124  return process
125 
126 # Prune gen particles with conditions applied in usual MiniAOD
128  if process.nanogenSequence.contains(process.mergedGenParticles):
129  raise ValueError("Applying the MiniAOD genParticle pruner to MiniAOD is redunant. " \
130  "Use a different customization.")
131  from PhysicsTools.PatAlgos.slimming.prunedGenParticles_cfi import prunedGenParticles
132  process.prunedGenParticles = prunedGenParticles.clone()
133  process.genParticleTable.src = "prunedGenParticles"
134  process.patJetPartons.particles = "prunedGenParticles"
135  process.nanogenSequence.insert(0, process.prunedGenParticles)
136  return process
137 
138 def setGenFullPrecision(process):
139  setGenPtPrecision(process, 23)
140  setGenEtaPrecision(process, 23)
141  setGenPhiPrecision(process, 23)
142 
143 def setGenPtPrecision(process, precision):
144  process.genParticleTable.variables.pt.precision = precision
145  process.genJetTable.variables.pt.precision = precision
146  process.metMCTable.variables.pt.precision = precision
147  return process
148 
149 def setGenEtaPrecision(process, precision):
150  process.genParticleTable.variables.eta.precision = precision
151  process.genJetTable.variables.eta.precision = precision
152  return process
153 
154 def setGenPhiPrecision(process, precision):
155  process.genParticleTable.variables.phi.precision = precision
156  process.genJetTable.variables.phi.precision = precision
157  process.metMCTable.variables.phi.precision = precision
158  return process
159 
160 def setLHEFullPrecision(process):
161  process.lheInfoTable.precision = 23
162  return process
163 
165  process.genWeightsTable.lheWeightPrecision = 23
166  return process
def pruneGenParticlesMini(process)
Definition: nanogen_cff.py:127
def customizeNanoGEN(process)
Definition: nanogen_cff.py:97
def setLHEFullPrecision(process)
Definition: nanogen_cff.py:160
def setGenPhiPrecision(process, precision)
Definition: nanogen_cff.py:154
def setGenPtPrecision(process, precision)
Definition: nanogen_cff.py:143
def nanoGenCommonCustomize(process)
Definition: nanogen_cff.py:44
def Var(expr, valtype, compression=None, doc=None, mcOnly=False, precision=-1)
Definition: common_cff.py:20
def pruneGenParticlesNano(process)
Definition: nanogen_cff.py:119
def setGenFullPrecision(process)
Definition: nanogen_cff.py:138
def setGenWeightsFullPrecision(process)
Definition: nanogen_cff.py:164
def customizeNanoGENFromMini(process)
Definition: nanogen_cff.py:62
def setGenEtaPrecision(process, precision)
Definition: nanogen_cff.py:149