1 |
peiffer |
1.1 |
# -*- coding: utf-8 -*-
|
2 |
|
|
import FWCore.ParameterSet.Config as cms
|
3 |
|
|
import os
|
4 |
|
|
|
5 |
|
|
process = cms.Process("MYPATTEST")
|
6 |
|
|
|
7 |
|
|
isttbar = False
|
8 |
|
|
issingletop = False
|
9 |
|
|
isvjets = False
|
10 |
|
|
isdata = True
|
11 |
|
|
nickname = "MC_TTbar"
|
12 |
|
|
#nickname= "DATA"
|
13 |
|
|
if nickname.startswith('MC_'): isdata = False
|
14 |
|
|
if nickname.find("_TTbar") > -1 or nickname.find('_ZP') > -1 : isttbar = True
|
15 |
|
|
if nickname.find("_ST_tchan") > -1 : issingletop = True
|
16 |
|
|
if nickname.find("_Wjets") > -1 or nickname.find('_Vqq') > -1 or nickname.find('_Wc') > -1 or nickname.find('_Zjets') > -1 : isvjets = True
|
17 |
|
|
tmp = os.getenv('CMSSW_BASE')
|
18 |
|
|
if 'CMSSW_4_1_' in tmp: cmssw_version = '41'
|
19 |
|
|
elif 'CMSSW_4_2_' in tmp: cmssw_version = '42'
|
20 |
|
|
elif 'CMSSW_4_4_' in tmp: cmssw_version = '44'
|
21 |
|
|
else: raise RuntimeError, 'could not determine cmssw version (cmssw correctly set up?)'
|
22 |
|
|
|
23 |
|
|
print "detected CMSSW version %s" % cmssw_version
|
24 |
|
|
if isdata: print "assuming to run over data"
|
25 |
|
|
else: print "assuming to run over MC; (single top, ttbar, vjets) = (", issingletop, isttbar, isvjets, ")"
|
26 |
|
|
|
27 |
|
|
# for global tag, see https://twiki.cern.ch/twiki/bin/view/CMSPublic/SWGuideFrontierConditions
|
28 |
|
|
def globaltag():
|
29 |
|
|
if isdata and cmssw_version=='41': return 'GR_R_41_V0::All'
|
30 |
|
|
if not isdata and cmssw_version=='41': return 'START41_V0::All'
|
31 |
|
|
if isdata and cmssw_version=='42': return 'GR_R_42_V25::All'
|
32 |
|
|
if not isdata and cmssw_version=='42': return 'START42_V17::All'
|
33 |
|
|
if isdata and cmssw_version=='44': return 'GR_R_44_V15::All'
|
34 |
|
|
if not isdata and cmssw_version=='44': return 'START44_V13::All'
|
35 |
|
|
|
36 |
|
|
raise RuntimeError, "could not determine global tag to use!"
|
37 |
|
|
|
38 |
|
|
|
39 |
|
|
process.load("FWCore.MessageLogger.MessageLogger_cfi")
|
40 |
|
|
process.MessageLogger.cerr.threshold = 'WARNING'
|
41 |
|
|
process.options = cms.untracked.PSet(
|
42 |
|
|
wantSummary = cms.untracked.bool(True)
|
43 |
|
|
)
|
44 |
|
|
|
45 |
|
|
process.source = cms.Source("PoolSource",
|
46 |
peiffer |
1.5 |
# fileNames = cms.untracked.vstring('dcap://dcache-cms-dcap.desy.de//pnfs/desy.de/cms/tier2/store/data/Run2011B/Jet/AOD/19Nov2011-v1/0000/985E38E8-3B15-E111-87A8-0015178C4970.root',
|
47 |
|
|
# 'dcap://dcache-cms-dcap.desy.de//pnfs/desy.de/cms/tier2/store/data/Run2011B/Jet/AOD/19Nov2011-v1/0000/4449A169-1A15-E111-B25F-0015178C4CE8.root'),
|
48 |
peiffer |
1.1 |
fileNames = cms.untracked.vstring('dcap://dcache-cms-dcap.desy.de//pnfs/desy.de/cms/tier2/store/mc/Fall11/TTJets_TuneZ2_7TeV-madgraph-tauola/AODSIM/PU_S6-START44_V5-v1/0000/68423CA2-BC04-E111-B90B-0018F3D0968C.root'),
|
49 |
peiffer |
1.5 |
# fileNames = cms.untracked.vstring('dcap://dcache-cms-dcap.desy.de//pnfs/desy.de/cms/tier2/store/mc/Fall11/TTJets_TuneZ2_7TeV-madgraph-tauola/AODSIM/PU_S6-START44_V5-v1/0000/7C7C20EE-BB04-E111-A833-001A92971BDA.root'),
|
50 |
peiffer |
1.1 |
#fileNames = cms.untracked.vstring('dcap://dcache-cms-dcap.desy.de//pnfs/desy.de/cms/tier2/store/mc/Fall11/W2Jets_TuneZ2_7TeV-madgraph-tauola/AODSIM/PU_S6_START44_V9B-v1/0000/9A4D5BEE-F536-E111-9F71-003048C692CA.root'),
|
51 |
|
|
#fileNames = cms.untracked.vstring('dcap://dcache-cms-dcap.desy.de//pnfs/desy.de/cms/tier2/store/mc/Fall11/QCD_Pt-80to170_EMEnriched_TuneZ2_7TeV-pythia6/AODSIM/PU_S6_START42_V14B-v2/0000/50404FDF-6214-E111-8248-78E7D164BFC8.root'),
|
52 |
|
|
skipEvents = cms.untracked.uint32(0),
|
53 |
|
|
duplicateCheckMode = cms.untracked.string('noDuplicateCheck')
|
54 |
|
|
)
|
55 |
|
|
|
56 |
|
|
|
57 |
|
|
process.maxEvents = cms.untracked.PSet( input = cms.untracked.int32(1000))
|
58 |
|
|
|
59 |
|
|
## Output Module Configuration (expects a path 'p')
|
60 |
|
|
from PhysicsTools.PatAlgos.patEventContent_cff import patEventContent
|
61 |
|
|
process.out = cms.OutputModule("PoolOutputModule",
|
62 |
|
|
fileName = cms.untracked.string('patTuple.root'),
|
63 |
|
|
# save only events passing the full path
|
64 |
|
|
SelectEvents = cms.untracked.PSet( SelectEvents = cms.vstring('p') ),
|
65 |
|
|
# save PAT Layer 1 output; you need a '*' to
|
66 |
|
|
# unpack the list of commands 'patEventContent'
|
67 |
|
|
outputCommands = cms.untracked.vstring('drop *' , *patEventContent)
|
68 |
|
|
)
|
69 |
|
|
|
70 |
|
|
if isdata and cmssw_version == '41': # note: there are no residual corrections in 42(?!);
|
71 |
|
|
# needs to be done with recorrector in SKITA
|
72 |
|
|
JECToProcess = cms.vstring(['L1FastJet','L2Relative', 'L3Absolute', 'L2L3Residual'])
|
73 |
|
|
else:
|
74 |
|
|
JECToProcess = cms.vstring(['L1FastJet','L2Relative', 'L3Absolute'])
|
75 |
|
|
|
76 |
|
|
# for some jet types, no residual correction is available; use l2l3 only:
|
77 |
|
|
#l2l3jec = cms.vstring(['L2Relative', 'L3Absolute'])
|
78 |
|
|
|
79 |
|
|
inputJetCorrLabel = ('AK5PFchs', JECToProcess)
|
80 |
|
|
|
81 |
|
|
# 1. general setup
|
82 |
|
|
process.load("Configuration.StandardSequences.Geometry_cff")
|
83 |
|
|
process.load("Configuration.StandardSequences.FrontierConditions_GlobalTag_cff")
|
84 |
|
|
process.GlobalTag.globaltag = globaltag()
|
85 |
|
|
print "using globaltag %s" % str(process.GlobalTag.globaltag)
|
86 |
|
|
process.load("Configuration.StandardSequences.MagneticField_cff")
|
87 |
|
|
|
88 |
|
|
|
89 |
|
|
|
90 |
|
|
# 3. configure PF2PAT
|
91 |
|
|
#process.load("PhysicsTools.PFCandProducer.PF2PAT_cff")
|
92 |
|
|
|
93 |
|
|
from RecoJets.JetProducers.kt4PFJets_cfi import kt4PFJets
|
94 |
|
|
process.kt6PFJets = kt4PFJets.clone()
|
95 |
|
|
process.kt6PFJets.rParam = 0.6
|
96 |
|
|
process.kt6PFJets.doRhoFastjet = True
|
97 |
|
|
|
98 |
|
|
process.kt6PFJetsRestrictedRho = process.kt6PFJets.clone()
|
99 |
|
|
process.kt6PFJetsRestrictedRho.Rho_EtaMax = cms.double(1.9)
|
100 |
|
|
process.kt6PFJetsRestrictedRho.Ghost_EtaMax = cms.double(2.5)
|
101 |
|
|
|
102 |
|
|
|
103 |
|
|
|
104 |
|
|
# Get a list of good primary vertices, in 42x, these are DAF vertices
|
105 |
|
|
from PhysicsTools.SelectorUtils.pvSelector_cfi import pvSelector
|
106 |
|
|
process.goodOfflinePrimaryVertices = cms.EDFilter(
|
107 |
|
|
"PrimaryVertexObjectFilter",
|
108 |
|
|
filterParams = pvSelector.clone( minNdof = cms.double(4.0), maxZ = cms.double(24.0) ),
|
109 |
|
|
src=cms.InputTag('offlinePrimaryVertices')
|
110 |
|
|
)
|
111 |
|
|
|
112 |
|
|
|
113 |
|
|
#top jets
|
114 |
|
|
|
115 |
peiffer |
1.2 |
process.load("UHHAnalysis.NtupleWriter.Topjets_cfi")
|
116 |
peiffer |
1.1 |
|
117 |
|
|
process.load("RecoJets.Configuration.GenJetParticles_cff")
|
118 |
|
|
from RecoJets.JetProducers.ca4GenJets_cfi import ca4GenJets
|
119 |
|
|
process.ca8GenJets = ca4GenJets.clone( rParam = cms.double(0.8) )
|
120 |
|
|
|
121 |
|
|
|
122 |
|
|
|
123 |
|
|
process.load("PhysicsTools.PatAlgos.patSequences_cff")
|
124 |
|
|
from PhysicsTools.PatAlgos.tools.trigTools import *
|
125 |
|
|
switchOnTrigger( process )
|
126 |
|
|
from PhysicsTools.PatAlgos.tools.pfTools import *
|
127 |
|
|
usePF2PAT(process, runPF2PAT=True, jetAlgo="AK5", runOnMC = not isdata)
|
128 |
|
|
process.pfPileUp.Enable = True
|
129 |
|
|
process.pfPileUp.Vertices = 'goodOfflinePrimaryVertices'
|
130 |
|
|
process.pfJets.doAreaFastjet = True
|
131 |
|
|
process.pfJets.doRhoFastjet = False
|
132 |
|
|
process.kt6PFJets.src = process.pfJets.src
|
133 |
|
|
process.kt6PFJetsRestrictedRho.src = process.pfJets.src
|
134 |
|
|
process.patJetCorrFactors.payload = 'AK5PFchs'
|
135 |
|
|
process.patJetCorrFactors.levels = JECToProcess
|
136 |
|
|
process.patJetCorrFactors.rho = cms.InputTag("kt6PFJets", "rho")
|
137 |
|
|
process.pfNoTau.enable = cms.bool(True) ##change this????
|
138 |
|
|
|
139 |
|
|
|
140 |
|
|
from PhysicsTools.PatAlgos.producersLayer1.photonProducer_cff import *
|
141 |
|
|
from PhysicsTools.PatAlgos.selectionLayer1.photonSelector_cfi import *
|
142 |
|
|
from PhysicsTools.PatAlgos.selectionLayer1.photonCountFilter_cfi import *
|
143 |
|
|
if isdata:
|
144 |
|
|
process.patDefaultSequence *= patPhotons
|
145 |
|
|
else: process.patDefaultSequence *= makePatPhotons
|
146 |
|
|
process.patDefaultSequence *= selectedPatPhotons
|
147 |
|
|
process.patDefaultSequence *= countPatPhotons
|
148 |
|
|
|
149 |
|
|
process.pfIsolatedMuons.isolationCut = cms.double(0.2)
|
150 |
|
|
|
151 |
|
|
|
152 |
|
|
process.patMuons.usePV = cms.bool(False)
|
153 |
|
|
process.patElectrons.usePV = cms.bool(False)
|
154 |
|
|
process.patJets.addTagInfos = cms.bool(True)
|
155 |
|
|
|
156 |
|
|
process.pfIsolatedMuonsLoose = process.pfIsolatedMuons.clone(isolationCut = cms.double(float("inf")))
|
157 |
|
|
process.patMuonsLoose = process.patMuons.clone(pfMuonSource = cms.InputTag("pfIsolatedMuonsLoose"), genParticleMatch = cms.InputTag("muonMatchLoose"))
|
158 |
|
|
# use pf isolation, but do not change matching (hell this "PAT" is really f*** up).
|
159 |
|
|
tmp = process.muonMatch.src
|
160 |
|
|
adaptPFMuons(process, process.patMuonsLoose, "")
|
161 |
|
|
process.muonMatch.src = tmp
|
162 |
|
|
|
163 |
|
|
|
164 |
|
|
process.pfSelectedElectrons.cut = cms.string('et > 15. & abs(eta) < 2.5')
|
165 |
|
|
process.pfSelectedMuons.cut = cms.string('pt > 10. & abs(eta) < 2.5')
|
166 |
|
|
process.pfSelectedPhotons.cut = cms.string('pt > 10.')
|
167 |
|
|
|
168 |
|
|
#jet and tau cuts are not propagated to ntuplewriter
|
169 |
|
|
process.selectedPatJets.cut = cms.string('pt > 10. & abs(eta) < 5.0')
|
170 |
|
|
process.selectedPatTaus.cut = cms.string('pt > 10. & abs(eta) < 5.0')
|
171 |
|
|
|
172 |
|
|
process.muonMatchLoose = process.muonMatch.clone(src = cms.InputTag("pfIsolatedMuonsLoose"))
|
173 |
|
|
|
174 |
|
|
process.pfIsolatedElectronsLoose = process.pfIsolatedElectrons.clone(isolationCut = cms.double(float("inf")))
|
175 |
|
|
process.patElectronsLoose = process.patElectrons.clone(pfElectronSource = cms.InputTag("pfIsolatedElectronsLoose"))
|
176 |
|
|
adaptPFElectrons(process, process.patElectronsLoose, "")
|
177 |
|
|
|
178 |
|
|
process.looseLeptonSequence = cms.Sequence(
|
179 |
|
|
process.pfIsolatedMuonsLoose +
|
180 |
|
|
process.muonMatchLoose +
|
181 |
|
|
process.patMuonsLoose +
|
182 |
|
|
process.pfIsolatedElectronsLoose +
|
183 |
|
|
process.patElectronsLoose
|
184 |
|
|
)
|
185 |
|
|
|
186 |
|
|
if isdata: process.looseLeptonSequence.remove(process.muonMatchLoose)
|
187 |
|
|
|
188 |
|
|
# note: we cannot use this for the moment: some modules require a primary vertex (even if is is fake!!)
|
189 |
|
|
#for module in (process.patMuons, process.patMuonsLoose, process.patElectrons, process.patElectronsLoose):
|
190 |
|
|
# module.pvSrc = "goodOfflinePrimaryVertices"
|
191 |
|
|
#for module in (process.impactParameterTagInfos, process.impactParameterTagInfosAOD, process.softMuonTagInfosAOD, process.softMuonTagInfos):
|
192 |
|
|
# module.primaryVertex = "goodOfflinePrimaryVertices"
|
193 |
|
|
#for module in (process.patJetCorrFactors,):
|
194 |
|
|
# module.primaryVertices = "goodOfflinePrimaryVertices"
|
195 |
|
|
#for module in (process.pfElectronsFromVertex,process.pfMuonsFromVertex):
|
196 |
|
|
# module.vertices = "goodOfflinePrimaryVertices"
|
197 |
|
|
|
198 |
|
|
|
199 |
|
|
|
200 |
|
|
process.load('CommonTools/RecoAlgos/HBHENoiseFilterResultProducer_cfi')
|
201 |
|
|
|
202 |
|
|
process.HBHENoiseFilterResultProducer.minRatio = cms.double(-999)
|
203 |
|
|
process.HBHENoiseFilterResultProducer.maxRatio = cms.double(999)
|
204 |
|
|
process.HBHENoiseFilterResultProducer.minHPDHits = cms.int32(17)
|
205 |
|
|
process.HBHENoiseFilterResultProducer.minRBXHits = cms.int32(999)
|
206 |
|
|
process.HBHENoiseFilterResultProducer.minHPDNoOtherHits = cms.int32(10)
|
207 |
|
|
process.HBHENoiseFilterResultProducer.minZeros = cms.int32(10)
|
208 |
|
|
process.HBHENoiseFilterResultProducer.minHighEHitTime = cms.double(-9999.0)
|
209 |
|
|
process.HBHENoiseFilterResultProducer.maxHighEHitTime = cms.double(9999.0)
|
210 |
|
|
process.HBHENoiseFilterResultProducer.maxRBXEMF = cms.double(-999.0)
|
211 |
|
|
process.HBHENoiseFilterResultProducer.minNumIsolatedNoiseChannels = cms.int32(999999)
|
212 |
|
|
process.HBHENoiseFilterResultProducer.minIsolatedNoiseSumE = cms.double(999999.)
|
213 |
|
|
process.HBHENoiseFilterResultProducer.minIsolatedNoiseSumEt = cms.double(999999.)
|
214 |
|
|
process.HBHENoiseFilterResultProducer.useTS4TS5 = cms.bool(True)
|
215 |
|
|
|
216 |
|
|
|
217 |
|
|
# add CA PF jets
|
218 |
|
|
process.caTopPFJets.src = process.pfJets.src
|
219 |
|
|
addJetCollection(process,
|
220 |
|
|
cms.InputTag('caTopPFJets'),
|
221 |
|
|
algoLabel = 'TopTag', typeLabel = 'PF', doJTA=False,
|
222 |
|
|
doBTagging=False, jetCorrLabel = inputJetCorrLabel, doType1MET=False,#use ak5 PFjet corrections as described in AN2011_194
|
223 |
|
|
genJetCollection = cms.InputTag("ca8GenJets"), doJetID = False)
|
224 |
|
|
|
225 |
|
|
# *** jet pruning ala Jott:
|
226 |
|
|
process.caPrunedPFJets.src = process.pfJets.src
|
227 |
|
|
addJetCollection(process,
|
228 |
|
|
cms.InputTag('caPrunedPFJets'),
|
229 |
|
|
algoLabel = 'Pruned', typeLabel = 'PF', doJTA=True,
|
230 |
|
|
doBTagging=True, jetCorrLabel = inputJetCorrLabel, doType1MET=False,#use ak5 PFjet corrections as described in AN2011_194
|
231 |
|
|
genJetCollection = cms.InputTag("ca8GenJets"), doJetID = False)
|
232 |
|
|
|
233 |
|
|
#configure jet parton flavour determination
|
234 |
|
|
for jetcollname in ('TopTagPF', 'PrunedPF'):
|
235 |
|
|
jetcoll = getattr(process, 'patJets'+jetcollname)
|
236 |
|
|
jetcoll.embedGenJetMatch = cms.bool(False)
|
237 |
|
|
jetcoll.addDiscriminators = cms.bool(False)
|
238 |
|
|
# jet flavor stuff: match also the top quark.
|
239 |
|
|
getattr(process, 'patJetPartonAssociation'+jetcollname).coneSizeToAssociate = cms.double(0.8)
|
240 |
|
|
getattr(process, 'patJetPartonAssociation'+jetcollname).doPriority = cms.bool(True)
|
241 |
|
|
getattr(process, 'patJetPartonAssociation'+jetcollname).priorityList = cms.vint32(6)
|
242 |
|
|
#4 = match heaviest flavor
|
243 |
|
|
getattr(process, 'patJetFlavourAssociation'+jetcollname).definition = cms.int32(4)
|
244 |
|
|
getattr(process, 'patJetFlavourAssociation'+jetcollname).physicsDefinition = cms.bool(False)
|
245 |
|
|
getattr(process, 'patJetPartonMatch'+jetcollname).mcPdgId = cms.vint32(1,2,3,4,5,6,21)
|
246 |
|
|
getattr(process, 'patJetPartonMatch'+jetcollname).maxDeltaR = cms.double(0.8)
|
247 |
|
|
|
248 |
|
|
for icorr in [process.patJetCorrFactorsTopTagPF,
|
249 |
|
|
process.patJetCorrFactorsPrunedPF
|
250 |
|
|
] :
|
251 |
|
|
icorr.rho = cms.InputTag("kt6PFJets", "rho")
|
252 |
|
|
|
253 |
|
|
# 5. compile main sequence
|
254 |
|
|
|
255 |
|
|
#process.patPF2PATSequence.replace(process.pfJets, process.kt6PFJets * process.kt6PFJetsRestrictedRho * process.pfJets)
|
256 |
|
|
process.patPF2PATSequence.replace(process.pfJets, process.kt6PFJets * process.kt6PFJetsRestrictedRho * process.pfJets * process.topjet_seq)
|
257 |
|
|
|
258 |
|
|
process.main_seq = cms.Sequence(process.goodOfflinePrimaryVertices)
|
259 |
|
|
|
260 |
|
|
process.main_seq *= process.patPF2PATSequence* process.looseLeptonSequence
|
261 |
|
|
|
262 |
|
|
process.main_genseq = cms.Sequence(process.genParticlesForJets * process.ca8GenJets * process.topjet_genseq )
|
263 |
|
|
#process.main_genseq = cms.Sequence(process.genParticlesForJets * process.ca8GenJets * process.subjet_genseq * process.topjet_genseq)
|
264 |
|
|
|
265 |
|
|
|
266 |
|
|
process.p = cms.Path()
|
267 |
|
|
|
268 |
|
|
if isdata: process.p *= process.HBHENoiseFilterResultProducer
|
269 |
|
|
if not isdata: process.p *= process.main_genseq
|
270 |
|
|
process.p *= process.main_seq
|
271 |
|
|
|
272 |
|
|
# Apply jet ID to all of the jets upstream. We aren't going to screw around
|
273 |
|
|
# with this, most likely. So, we don't really to waste time with it
|
274 |
|
|
# at the analysis level.
|
275 |
|
|
from PhysicsTools.SelectorUtils.pfJetIDSelector_cfi import pfJetIDSelector
|
276 |
|
|
process.goodPatJets = cms.EDFilter("PFJetIDSelectionFunctorFilter",
|
277 |
|
|
filterParams = pfJetIDSelector.clone(),
|
278 |
|
|
src = cms.InputTag("patJets")
|
279 |
|
|
)
|
280 |
|
|
|
281 |
|
|
process.goodPatJetsPrunedPF = cms.EDFilter("PFJetIDSelectionFunctorFilter",
|
282 |
|
|
filterParams = pfJetIDSelector.clone(),
|
283 |
|
|
src = cms.InputTag("patJetsPrunedPF")
|
284 |
|
|
)
|
285 |
|
|
|
286 |
|
|
|
287 |
|
|
process.goodPatJetsTopTagPF = cms.EDFilter("PFJetIDSelectionFunctorFilter",
|
288 |
|
|
filterParams = pfJetIDSelector.clone(),
|
289 |
|
|
src = cms.InputTag("patJetsTopTagPF")
|
290 |
|
|
)
|
291 |
|
|
|
292 |
|
|
|
293 |
|
|
|
294 |
|
|
#NtupleWriter
|
295 |
|
|
|
296 |
|
|
process.MyNtuple = cms.EDAnalyzer('NtupleWriter',
|
297 |
|
|
fileName = cms.string('ttbar.root'),
|
298 |
|
|
doElectrons = cms.bool(True),
|
299 |
|
|
doMuons = cms.bool(True),
|
300 |
|
|
doTaus = cms.bool(True),
|
301 |
|
|
doJets = cms.bool(True),
|
302 |
|
|
doTopJets = cms.bool(True),
|
303 |
|
|
doPhotons = cms.bool(True),
|
304 |
|
|
doMET = cms.bool(True),
|
305 |
|
|
doPV = cms.bool(True),
|
306 |
|
|
doGenInfo = cms.bool(not isdata),
|
307 |
peiffer |
1.5 |
doTrigger = cms.bool(True),
|
308 |
peiffer |
1.1 |
electron_sources = cms.vstring("patElectrons","patElectronsLoose"),
|
309 |
|
|
muon_sources = cms.vstring("patMuons","patMuonsLoose"),
|
310 |
|
|
tau_sources = cms.vstring("patTaus"),
|
311 |
|
|
tau_ptmin = cms.double(10.0),
|
312 |
|
|
tau_etamax = cms.double(5.0),
|
313 |
|
|
jet_sources = cms.vstring("goodPatJets"),
|
314 |
|
|
jet_ptmin = cms.double(10.0),
|
315 |
|
|
jet_etamax = cms.double(5.0),
|
316 |
|
|
topjet_sources = cms.vstring("goodPatJetsTopTagPF","goodPatJetsPrunedPF"),
|
317 |
|
|
topjet_ptmin = cms.double(150.0),
|
318 |
|
|
topjet_etamax = cms.double(5.0),
|
319 |
|
|
photon_sources = cms.vstring("patPhotons"),
|
320 |
|
|
met_sources = cms.vstring("patMETs"),
|
321 |
|
|
pv_sources = cms.vstring("goodOfflinePrimaryVertices", "offlinePrimaryVertices"),
|
322 |
|
|
trigger_prefixes = cms.vstring(#"HLT_IsoMu", "HLT_Mu",
|
323 |
|
|
#"HLT_L1SingleMu", "HLT_L2Mu",
|
324 |
|
|
#"HLT_Ele",
|
325 |
|
|
"HLT_Jet", "HLT_HT",
|
326 |
|
|
#"HLT_DoubleMu", "HLT_DoubleEle"
|
327 |
|
|
),
|
328 |
|
|
|
329 |
|
|
)
|
330 |
|
|
|
331 |
|
|
#process.TFileService = cms.Service("TFileService",
|
332 |
|
|
# fileName = cms.string('tree3.root')
|
333 |
|
|
#)
|
334 |
|
|
process.p *=process.goodPatJets
|
335 |
|
|
process.p *=process.goodPatJetsPrunedPF
|
336 |
|
|
process.p *=process.goodPatJetsTopTagPF
|
337 |
|
|
process.p *= process.MyNtuple
|
338 |
|
|
|
339 |
|
|
#comment this line in to produce pattuples:
|
340 |
|
|
#process.outpath = cms.EndPath(process.out)
|
341 |
|
|
|
342 |
|
|
process.options.wantSummary = True
|
343 |
|
|
|
344 |
|
|
#open('dump.py', 'w').write(process.dumpPython())
|