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