CMS 3D CMS Logo

btvMC_cff.py
Go to the documentation of this file.
1 import FWCore.ParameterSet.Config as cms
2 
4 from PhysicsTools.NanoAOD.simpleCandidateFlatTableProducer_cfi import simpleCandidateFlatTableProducer
5 from PhysicsTools.NanoAOD.jetsAK8_cff import fatJetTable, subJetTable
6 from PhysicsTools.NanoAOD.jetsAK4_Puppi_cff import jetPuppiTable
8 # when storing the flat table, always do "T"able in the naming convention
9 
10 
11 finalGenParticles.select +=[
12  "keep (4 <= abs(pdgId) <= 5) && statusFlags().isLastCopy()", # BTV: keep b/c quarks in their last copy
13  "keep (abs(pdgId) == 310 || abs(pdgId) == 3122) && statusFlags().isLastCopy()", # BTV: keep K0s and Lambdas in their last copy
14  ]
15 
16 btvGenTable = cms.EDProducer(
17  "SimpleGenParticleFlatTableProducer",
18  src=finalGenParticles.src,
19  name= cms.string("GenPart"),
20  doc = cms.string("interesting gen particles "),
21  singleton=cms.bool(False),
22  variables=
23  cms.PSet(
24  genParticleTable.variables,
25  vx = Var("vx", "float", doc="x coordinate of vertex position"),
26  vy = Var("vy", "float", doc="y coordinate of vertex position"),
27  vz = Var("vz", "float", doc="z coordinate of vertex position"),
28  genPartIdxMother2 = Var("?numberOfMothers>1?motherRef(1).key():-1", "int", doc="index of the second mother particle, if valid")
29  ))
30 genParticleTablesTask.replace(genParticleTable,btvGenTable)
31 btvMCTable = cms.EDProducer("BTVMCFlavourTableProducer",name=jetPuppiTable.name,src=cms.InputTag("linkedObjects","jets"),genparticles=cms.InputTag("prunedGenParticles"))
32 
33 btvAK4JetExtTable = cms.EDProducer(
34  "SimplePATJetFlatTableProducer",
35  src=jetPuppiTable.src,
36  cut=jetPuppiTable.cut,
37  name=jetPuppiTable.name,
38  doc=jetPuppiTable.doc,
39  singleton=cms.bool(False), # the number of entries is variable
40  extension=cms.bool(True), # this is the extension table for Jets
41  variables=cms.PSet(
42  nBHadrons=Var("jetFlavourInfo().getbHadrons().size()",
43  int,
44  doc="number of b-hadrons"),
45  nCHadrons=Var("jetFlavourInfo().getcHadrons().size()",
46  int,
47  doc="number of c-hadrons"),
48  ))
49 
50 btvSubJetMCExtTable = cms.EDProducer(
51  "SimplePATJetFlatTableProducer",
52  src = subJetTable.src,
53  cut = subJetTable.cut,
54  name = subJetTable.name,
55  doc=subJetTable.doc,
56  singleton = cms.bool(False),
57  extension = cms.bool(True),
58  variables = cms.PSet(
59  subGenJetAK8Idx = Var("?genJetFwdRef().backRef().isNonnull()?genJetFwdRef().backRef().key():-1",
60  int,
61  doc="index of matched gen Sub jet"),
62  )
63  )
64 genJetsAK8Constituents = cms.EDProducer("GenJetPackedConstituentPtrSelector",
65  src = cms.InputTag("slimmedGenJetsAK8"),
66  cut = cms.string("pt > 100.")
67  )
68 
69 
70 genJetsAK4Constituents = cms.EDProducer("GenJetPackedConstituentPtrSelector",
71  src = cms.InputTag("slimmedGenJets"),
72  cut = cms.string("pt > 20")
73  )
74 
75 
76 
77 
78 
79 ak4onlygenJetsConstituents = cms.EDProducer("PackedGenParticlePtrMerger", src = cms.VInputTag(cms.InputTag("genJetsAK4Constituents", "constituents")), skipNulls = cms.bool(True), warnOnSkip = cms.bool(True))
80 
81 ak8onlygenJetsConstituents = cms.EDProducer("PackedGenParticlePtrMerger", src = cms.VInputTag(cms.InputTag("genJetsAK8Constituents", "constituents")), skipNulls = cms.bool(True), warnOnSkip = cms.bool(True))
82 ak4ak8genJetsConstituents = cms.EDProducer("PackedGenParticlePtrMerger", src = cms.VInputTag(cms.InputTag("genJetsAK4Constituents", "constituents"), cms.InputTag("genJetsAK8Constituents", "constituents")), skipNulls = cms.bool(True), warnOnSkip = cms.bool(True))
83 
84 allPFParticleTable = cms.EDProducer("SimpleCandidateFlatTableProducer",
85  src = cms.InputTag("packedGenParticles"),
86  cut = cms.string(""), #we should not filter after pruning
87  name= cms.string("GenCands"),
88  doc = cms.string("interesting gen particles from PF candidates"),
89  singleton = cms.bool(False), # the number of entries is variable
90  extension = cms.bool(False), # this is the main table for the AK8 constituents
91  variables = cms.PSet(CandVars
92  )
93  )
94 ak4onlygenJetsParticleTable = cms.EDProducer("SimpleCandidateFlatTableProducer",
95  src = cms.InputTag("ak4onlygenJetsConstituents"),
96  cut = cms.string(""), #we should not filter after pruning
97  name= cms.string("GenCands"),
98  doc = cms.string("interesting gen particles from AK4 jets"),
99  singleton = cms.bool(False), # the number of entries is variable
100  extension = cms.bool(False), # this is the main table for the AK8 constituents
101  variables = cms.PSet(CandVars
102  )
103  )
104 ak8onlygenJetsParticleTable = cms.EDProducer("SimpleCandidateFlatTableProducer",
105  src = cms.InputTag("ak8onlygenJetsConstituents"),
106  cut = cms.string(""), #we should not filter after pruning
107  name= cms.string("GenCands"),
108  doc = cms.string("interesting gen particles from AK8 jets"),
109  singleton = cms.bool(False), # the number of entries is variable
110  extension = cms.bool(False), # this is the main table for the AK8 constituents
111  variables = cms.PSet(CandVars
112  )
113  )
114 ak4ak8genJetsParticleTable = cms.EDProducer("SimpleCandidateFlatTableProducer",
115  src = cms.InputTag("ak4ak8genJetsConstituents"),
116  cut = cms.string(""), #we should not filter after pruning
117  name= cms.string("GenCands"),
118  doc = cms.string("interesting gen particles from AK4, AK8 jets"),
119  singleton = cms.bool(False), # the number of entries is variable
120  extension = cms.bool(False), # this is the main table for the AK8 constituents
121  variables = cms.PSet(CandVars
122  )
123  )
124 ak8onlygenAK8ConstituentsTable = cms.EDProducer("GenJetConstituentTableProducer",
125  candidates = cms.InputTag("ak8onlygenJetsConstituents"),
126  jets = cms.InputTag("genJetsAK8Constituents"), # Note: The name has "Constituents" in it, but these are the jets
127  name = cms.string("GenFatJetCands"),
128  nameSV = cms.string("GenFatJetSVs"),
129  idx_name = cms.string("pFCandsIdx"),
130  idx_nameSV = cms.string("sVIdx"),
131  readBtag = cms.bool(False))
132 ak4onlygenAK4ConstituentsTable = cms.EDProducer("GenJetConstituentTableProducer",
133  candidates = cms.InputTag("ak4onlygenJetsConstituents"),
134  jets = cms.InputTag("genJetsAK4Constituents"), # Note: The name has "Constituents" in it, but these are the jets
135  name = cms.string("GenJetCands"),
136  nameSV = cms.string("GenJetSVs"),
137  idx_name = cms.string("pFCandsIdx"),
138  idx_nameSV = cms.string("sVIdx"),
139  readBtag = cms.bool(False))
140 ak4ak8genAK4ConstituentsTable = cms.EDProducer("GenJetConstituentTableProducer",
141  candidates = cms.InputTag("ak4ak8genJetsConstituents"),
142  jets = cms.InputTag("genJetsAK4Constituents"), # Note: The name has "Constituents" in it, but these are the jets
143  name = cms.string("GenJetCands"),
144  nameSV = cms.string("GenJetSVs"),
145  idx_name = cms.string("pFCandsIdx"),
146  idx_nameSV = cms.string("sVIdx"),
147  readBtag = cms.bool(False))
148 
149 ak4ak8genAK8ConstituentsTable = cms.EDProducer("GenJetConstituentTableProducer",
150  candidates = cms.InputTag("ak4ak8genJetsConstituents"),
151  jets = cms.InputTag("genJetsAK8Constituents"), # Note: The name has "Constituents" in it, but these are the jets
152  name = cms.string("GenFatJetCands"),
153  nameSV = cms.string("GenFatJetSVs"),
154  idx_name = cms.string("pFCandsIdx"),
155  idx_nameSV = cms.string("sVIdx"),
156  readBtag = cms.bool(False))
157 btvAK4MCSequence = cms.Sequence(btvGenTable+btvAK4JetExtTable+btvMCTable)
158 btvAK8MCSequence = cms.Sequence(btvGenTable+btvSubJetMCExtTable)
159 #PF Cands of AK4 only , with cross linking to AK4 jets
160 ak4onlyPFCandsMCSequence=cms.Sequence(genJetsAK4Constituents+ak4onlygenJetsConstituents+ak4onlygenJetsParticleTable+ak4onlygenAK4ConstituentsTable)+btvAK4MCSequence
161 #PF Cands of AK8 only , with cross linking to AK8 jets
162 ak8onlyPFCandsMCSequence=cms.Sequence(genJetsAK8Constituents+ak8onlygenJetsConstituents+ak8onlygenJetsParticleTable+ak8onlygenAK8ConstituentsTable)+btvAK8MCSequence
163 #PF Cands of AK4+AK8, with cross linking to AK4,AK8 jets
164 ak4ak8PFCandsMCSequence=cms.Sequence(genJetsAK4Constituents+genJetsAK8Constituents+ak4ak8genJetsConstituents+ak4ak8genJetsParticleTable+ak4ak8genAK4ConstituentsTable+ak4ak8genAK8ConstituentsTable)+btvAK4MCSequence+btvAK8MCSequence
165 #All PFCands, with cross linking to AK4,AK8 jets
166 allPFPFCandsMCSequence=cms.Sequence(genJetsAK4Constituents+genJetsAK8Constituents+ak4ak8genJetsConstituents+allPFParticleTable+ak4ak8genAK4ConstituentsTable+ak4ak8genAK8ConstituentsTable)+btvAK4MCSequence+btvAK8MCSequence
167 
168 
169 
170 
def Var(expr, valtype, doc=None, precision=-1, lazyEval=False)
Definition: common_cff.py:17