CMS 3D CMS Logo

Public Member Functions | Private Attributes | Static Private Attributes

jetTools::AddJetCollection Class Reference

Inherits FWCore::GuiBrowsers::ConfigToolBase::ConfigToolBase.

List of all members.

Public Member Functions

def __call__
def __init__
def getDefaultParameters
def toolCode

Private Attributes

 _comment
 _parameters

Static Private Attributes

tuple _defaultParameters = dicttypes.SortedKeysDict()
string _label = 'addJetCollection'

Detailed Description

Add a new collection of jets. Takes the configuration from the
already configured standard jet collection as starting point;
replaces before calling addJetCollection will also affect the
new jet collections

Definition at line 228 of file jetTools.py.


Constructor & Destructor Documentation

def jetTools::AddJetCollection::__init__ (   self)

Definition at line 238 of file jetTools.py.

00239                       :
00240         ConfigToolBase.__init__(self)
00241         self.addParameter(self._defaultParameters,'jetCollection',self._defaultValue,'Input jet collection', cms.InputTag)
00242         self.addParameter(self._defaultParameters,'algoLabel',self._defaultValue, "label to indicate the jet algorithm (e.g.'AK5')",str)
00243         self.addParameter(self._defaultParameters,'typeLabel',self._defaultValue, "label to indicate the type of constituents (e.g. 'Calo', 'Pflow', 'Jpt', ...)",str)
00244         self.addParameter(self._defaultParameters,'btagInfo',['impactParameterTagInfos','secondaryVertexTagInfos','softMuonTagInfos','secondaryVertexNegativeTagInfos','secondaryVertexNegativeTagInfos','inclusiveSecondaryVertexFinderTagInfos','softElectronTagInfos'],"input btag info",allowedValues=['impactParameterTagInfos','secondaryVertexTagInfos','softMuonTagInfos','secondaryVertexNegativeTagInfos','inclusiveSecondaryVertexFinderTagInfos','softElectronTagInfos'],Type=list)
00245         self.addParameter(self._defaultParameters,'btagdiscriminators',['jetBProbabilityBJetTags', 'jetProbabilityBJetTags','trackCountingHighPurBJetTags','trackCountingHighEffBJetTags','simpleSecondaryVertexHighEffBJetTags','simpleSecondaryVertexHighPurBJetTags','combinedSecondaryVertexBJetTags','combinedSecondaryVertexMVABJetTags','softMuonBJetTags','softMuonByPtBJetTags','softMuonByIP3dBJetTags','simpleSecondaryVertexNegativeHighEffBJetTags','simpleSecondaryVertexNegativeHighPurBJetTags','negativeTrackCountingHighEffJetTags','negativeTrackCountingHighPurJetTags','combinedInclusiveSecondaryVertexBJetTags','combinedMVABJetTags'],"input btag discriminators", allowedValues=['jetBProbabilityBJetTags', 'jetProbabilityBJetTags', 'trackCountingHighPurBJetTags', 'trackCountingHighEffBJetTags', 'simpleSecondaryVertexHighEffBJetTags','simpleSecondaryVertexHighPurBJetTags','combinedSecondaryVertexBJetTags','combinedSecondaryVertexMVABJetTags','softMuonBJetTags','softMuonByPtBJetTags','softMuonByIP3dBJetTags','simpleSecondaryVertexNegativeHighEffBJetTags','simpleSecondaryVertexNegativeHighPurBJetTags','negativeTrackCountingHighEffJetTags','negativeTrackCountingHighPurJetTags','combinedInclusiveSecondaryVertexBJetTags','combinedMVABJetTags'],Type=list)
00246         self.addParameter(self._defaultParameters,'doJTA',True, "run b tagging sequence for new jet collection and add it to the new pat jet collection")
00247         self.addParameter(self._defaultParameters,'doBTagging',True, 'run JetTracksAssociation and JetCharge and add it to the new pat jet collection (will autom. be true if doBTagging is set to true)')
00248         self.addParameter(self._defaultParameters,'jetCorrLabel',None, "payload and list of new jet correction labels, such as (\'AK5Calo\',[\'L2Relative\', \'L3Absolute\'])", tuple,acceptNoneValue=True )
00249         self.addParameter(self._defaultParameters,'doType1MET',True, "if jetCorrLabel is not 'None', set this to 'True' to redo the Type1 MET correction for the new jet colllection; at the moment it must be 'False' for non CaloJets otherwise the JetMET POG module crashes. ")
00250         self.addParameter(self._defaultParameters,'doL1Cleaning',True, "copy also the producer modules for cleanLayer1 will be set to 'True' automatically when doL1Counters is 'True'")
00251         self.addParameter(self._defaultParameters,'doL1Counters',False, "copy also the filter modules that accept/reject the event looking at the number of jets")
00252         self.addParameter(self._defaultParameters,'genJetCollection',cms.InputTag("ak5GenJets"), "GenJet collection to match to")
00253         self.addParameter(self._defaultParameters,'doJetID',True, "add jetId variables to the added jet collection?")
00254         self.addParameter(self._defaultParameters,'jetIdLabel',"ak5", " specify the label prefix of the xxxJetID object; in general it is the jet collection tag like ak5, kt4 sc5, aso. For more information have a look to SWGuidePATTools#add_JetCollection")
00255         self.addParameter(self._defaultParameters,'standardAlgo',"AK5", "standard algorithm label of the collection from which the clones for the new jet collection will be taken from (note that this jet collection has to be available in the event before hand)")
00256         self.addParameter(self._defaultParameters,'standardType',"Calo", "standard constituent type label of the collection from which the clones for the new jet collection will be taken from (note that this jet collection has to be available in the event before hand)")
00257         self.addParameter(self._defaultParameters, 'outputModules', ['out'], "output module labels, empty list of label indicates no output, default: ['out']")
00258 
00259         self._parameters=copy.deepcopy(self._defaultParameters)
00260         self._comment = ""


Member Function Documentation

def jetTools::AddJetCollection::__call__ (   self,
  process,
  jetCollection = None,
  algoLabel = None,
  typeLabel = None,
  doJTA = None,
  doBTagging = None,
  jetCorrLabel = None,
  doType1MET = None,
  doL1Cleaning = None,
  doL1Counters = None,
  genJetCollection = None,
  doJetID = None,
  jetIdLabel = None,
  outputModule = None,
  outputModules = None,
  btagInfo = None,
  btagdiscriminators = None 
)

Definition at line 264 of file jetTools.py.

00283                                        :
00284 
00285         ## stop processing if 'outputModule' exists and show the new alternative
00286         if  not outputModule is None:
00287             depricatedOptionOutputModule(self)
00288         if jetCollection  is None:
00289             jetCollection=self._defaultParameters['jetCollection'].value
00290         if algoLabel is None:
00291             algoLabel=self._defaultParameters['algoLabel'].value
00292         if typeLabel is None:
00293             typeLabel=self._defaultParameters['typeLabel'].value
00294         if doJTA is None:
00295             doJTA=self._defaultParameters['doJTA'].value
00296         if doBTagging is None:
00297             doBTagging=self._defaultParameters['doBTagging'].value
00298         if jetCorrLabel  is None:
00299             jetCorrLabel=self._defaultParameters['jetCorrLabel'].value
00300         if doType1MET  is None:
00301             doType1MET=self._defaultParameters['doType1MET'].value
00302         if doL1Cleaning is None:
00303             doL1Cleaning=self._defaultParameters['doL1Cleaning'].value
00304         if doL1Counters  is None:
00305             doL1Counters=self._defaultParameters['doL1Counters'].value
00306         if genJetCollection  is None:
00307             genJetCollection=self._defaultParameters['genJetCollection'].value
00308         if doJetID  is None:
00309             doJetID=self._defaultParameters['doJetID'].value
00310         if jetIdLabel  is None:
00311             jetIdLabel=self._defaultParameters['jetIdLabel'].value
00312         if outputModules is None:
00313             outputModules=self._defaultParameters['outputModules'].value
00314         if  btagInfo is None:
00315             btagInfo=self._defaultParameters['btagInfo'].value
00316         if  btagdiscriminators is None:
00317             btagdiscriminators=self._defaultParameters['btagdiscriminators'].value
00318 
00319         self.setParameter('jetCollection',jetCollection)
00320         self.setParameter('algoLabel',algoLabel)
00321         self.setParameter('typeLabel',typeLabel)
00322         self.setParameter('doJTA',doJTA)
00323         self.setParameter('doBTagging',doBTagging)
00324         self.setParameter('jetCorrLabel',jetCorrLabel)
00325         self.setParameter('doType1MET',doType1MET)
00326         self.setParameter('doL1Cleaning',doL1Cleaning)
00327         self.setParameter('doL1Counters',doL1Counters)
00328         self.setParameter('genJetCollection',genJetCollection)
00329         self.setParameter('doJetID',doJetID)
00330         self.setParameter('jetIdLabel',jetIdLabel)
00331         self.setParameter('outputModules',outputModules)
00332         self.setParameter('btagInfo',btagInfo)
00333         self.setParameter('btagdiscriminators',btagdiscriminators)
00334 
00335         self.apply(process)

def jetTools::AddJetCollection::getDefaultParameters (   self)

Definition at line 261 of file jetTools.py.

00262                                   :
00263         return self._defaultParameters

def jetTools::AddJetCollection::toolCode (   self,
  process 
)

Definition at line 336 of file jetTools.py.

00337                                :
00338 
00339         jetCollection=self._parameters['jetCollection'].value
00340         algoLabel=self._parameters['algoLabel'].value
00341         typeLabel=self._parameters['typeLabel'].value
00342         doJTA=self._parameters['doJTA'].value
00343         doBTagging=self._parameters['doBTagging'].value
00344         jetCorrLabel=self._parameters['jetCorrLabel'].value
00345         doType1MET =self._parameters['doType1MET'].value
00346         doL1Cleaning=self._parameters['doL1Cleaning'].value
00347         doL1Counters=self._parameters['doL1Counters'].value
00348         genJetCollection=self._parameters['genJetCollection'].value
00349         doJetID=self._parameters['doJetID'].value
00350         jetIdLabel=self._parameters['jetIdLabel'].value
00351         outputModules=self._parameters['outputModules'].value
00352         btagInfo=self._parameters['btagInfo'].value
00353         btagdiscriminators=self._parameters['btagdiscriminators'].value
00354 
00355 
00356         ## create old module label from standardAlgo
00357         ## and standardType and return
00358         def oldLabel(prefix=''):
00359             return jetCollectionString(prefix, '', '')
00360 
00361         ## create new module label from old module
00362         ## label and return
00363         def newLabel(oldLabel):
00364             newLabel=oldLabel
00365             oldLabel=oldLabel+algoLabel+typeLabel
00366             return oldLabel
00367 
00368         ## clone module and add it to the patDefaultSequence
00369         def addClone(hook, **replaceStatements):
00370             ## create a clone of the hook with corresponding
00371             ## parameter replacements
00372             newModule = getattr(process, hook).clone(**replaceStatements)
00373             ## add the module to the sequence
00374             addModuleToSequence(hook, newModule)
00375 
00376         ## add module to the patDefaultSequence
00377         def addModuleToSequence(hook, newModule):
00378             hookModule = getattr(process, hook)
00379             ## add the new module with standardAlgo &
00380             ## standardType replaced in module label
00381             setattr( process, newLabel(hook), newModule)
00382             ## add new module to default sequence
00383             ## just behind the hookModule
00384             process.patDefaultSequence.replace( hookModule, hookModule*newModule )
00385 
00386         ## add a clone of patJets
00387         addClone(oldLabel(), jetSource = jetCollection)
00388         ## add a clone of selectedPatJets
00389         addClone(oldLabel('selected'), src=cms.InputTag(newLabel(oldLabel())))
00390         ## add a clone of cleanPatJets
00391         if( doL1Cleaning ):
00392             addClone(oldLabel('clean'), src=cms.InputTag(newLabel(oldLabel('selected'))))
00393         ## add a clone of countPatJets
00394         if( doL1Counters ):
00395             if( doL1Cleaning ):
00396                 addClone(oldLabel('count'), src=cms.InputTag(newLabel(oldLabel('clean'))))
00397             else:
00398                 addClone(oldLabel('count'), src=cms.InputTag(newLabel(oldLabel('selected'))))
00399 
00400         ## get attributes of new module
00401         l1Jets = getattr(process, newLabel(oldLabel()))
00402 
00403         ## add a clone of gen jet matching
00404         addClone('patJetPartonMatch', src = jetCollection)
00405         addClone('patJetGenJetMatch', src = jetCollection, matched = genJetCollection)
00406 
00407         ## add a clone of parton and flavour associations
00408         addClone('patJetPartonAssociation', jets = jetCollection)
00409         addClone('patJetFlavourAssociation', srcByReference = cms.InputTag(newLabel('patJetPartonAssociation')))
00410 
00411         ## fix label for input tag
00412         def fixInputTag(x): x.setModuleLabel(newLabel(x.moduleLabel))
00413         ## fix label for vector of input tags
00414         def fixVInputTag(x): x[0].setModuleLabel(newLabel(x[0].moduleLabel))
00415 
00416         ## provide allLayer1Jet inputs with individual labels
00417         fixInputTag(l1Jets.genJetMatch)
00418         fixInputTag(l1Jets.genPartonMatch)
00419         fixInputTag(l1Jets.JetPartonMapSource)
00420 
00421         ## make VInputTag from strings
00422         def vit(*args) : return cms.VInputTag( *[ cms.InputTag(x) for x in args ] )
00423 
00424         if (doJTA or doBTagging):
00425             ## add clone of jet track association
00426             process.load("RecoJets.JetAssociationProducers.ak5JTA_cff")
00427             from RecoJets.JetAssociationProducers.ak5JTA_cff import ak5JetTracksAssociatorAtVertex
00428             ## add jet track association module to processes
00429             jtaLabel = 'jetTracksAssociatorAtVertex'+algoLabel+typeLabel
00430             setattr( process, jtaLabel, ak5JetTracksAssociatorAtVertex.clone(jets = jetCollection) )
00431             process.patDefaultSequence.replace(process.patJetCharge, getattr(process,jtaLabel)+process.patJetCharge)
00432             l1Jets.trackAssociationSource = cms.InputTag(jtaLabel)
00433             addClone('patJetCharge', src=cms.InputTag(jtaLabel)),
00434             fixInputTag(l1Jets.jetChargeSource)
00435         else:
00436             ## switch embedding of track association and jet
00437             ## charge estimate to 'False'
00438             l1Jets.addAssociatedTracks = False
00439             l1Jets.addJetCharge = False
00440 
00441         if (doBTagging):
00442             ## define postfixLabel
00443             postfixLabel=algoLabel+typeLabel
00444             ## add b tagging sequence
00445             (btagSeq, btagLabels) = runBTagging(process, jetCollection, postfixLabel,"", btagInfo,btagdiscriminators)
00446             ## add b tagging sequence before running the allLayer1Jets modules
00447             process.patDefaultSequence.replace(getattr(process,jtaLabel), getattr(process,jtaLabel)+btagSeq)
00448             ## replace corresponding tags for pat jet production
00449             l1Jets.trackAssociationSource = cms.InputTag(btagLabels['jta'])
00450             l1Jets.tagInfoSources = cms.VInputTag( *[ cms.InputTag(x) for x in btagLabels['tagInfos'] ] )
00451             l1Jets.discriminatorSources = cms.VInputTag( *[ cms.InputTag(x) for x in btagLabels['jetTags']  ] )
00452         else:
00453             ## switch general b tagging info switch off
00454             l1Jets.addBTagInfo = False
00455             ## adjust output
00456             if len(outputModules) > 0:
00457                 for outMod in outputModules:
00458                     if hasattr(process,outMod):
00459                         getattr(process,outMod).outputCommands.append("drop *_"+newLabel(oldLabel('selected'))+"_tagInfos_*")
00460                     else:
00461                         raise KeyError, "process has no OutModule named", outMod
00462 
00463         if (doJetID):
00464             l1Jets.addJetID = cms.bool(True)
00465             jetIdLabelNew = jetIdLabel + 'JetID'
00466             l1Jets.jetIDMap = cms.InputTag( jetIdLabelNew )
00467         else :
00468             l1Jets.addJetID = cms.bool(False)
00469 
00470         if (jetCorrLabel != None):
00471             ## add clone of jet energy corrections;
00472             ## catch a couple of exceptions first
00473             if (jetCorrLabel == False ):
00474                 raise ValueError, "In addJetCollection 'jetCorrLabel' must be set to 'None', not 'False'"
00475             if (jetCorrLabel == "None"):
00476                 raise ValueError, "In addJetCollection 'jetCorrLabel' must be set to 'None' (without quotes)"
00477             ## check for the correct format
00478             if type(jetCorrLabel) != type(('AK5Calo',['L2Relative'])):
00479                 raise ValueError, "In addJetCollection 'jetCorrLabel' must be 'None', or of type ('payload',['correction1', 'correction2'])"
00480 
00481             ## add clone of jetCorrFactors
00482             addClone('patJetCorrFactors', src = jetCollection)
00483             switchJetCorrLevels(process, jetCorrLabel = jetCorrLabel, postfix=algoLabel+typeLabel)
00484             getattr(process, newLabel('patJets')).jetCorrFactorsSource = cms.VInputTag(  cms.InputTag(newLabel('patJetCorrFactors')) )
00485 
00486             ## find out type of jet collection, switch type1MET corrections off for JPTJets
00487             jetCollType = ''
00488             if   ( 'CaloJets' in jetCollection.getModuleLabel() ):
00489                 jetCollType = 'Calo'
00490             elif ( 'PFJets' in jetCollection.getModuleLabel() or jetCollection.getModuleLabel().startswith('pfNo') or jetCollection.getModuleLabel() == 'particleFlow'):
00491                 jetCollType = 'PF'
00492             else:
00493                 print '============================================='
00494                 print 'Type1MET corrections are switched off for    '
00495                 print 'JPT Jets. Users are recommened to use tcMET  '
00496                 print 'together with JPT jets.                      '
00497                 print '============================================='
00498                 doType1MET=False
00499 
00500             ## add a clone of the type1MET correction for the new jet collection
00501             if (doType1MET):
00502                 ## create jet correctors for MET corrections
00503                 from JetMETCorrections.Configuration.JetCorrectionServicesAllAlgos_cff import ak5PFL1Fastjet, ak5PFL1Offset, ak5PFL2Relative, ak5PFL3Absolute, ak5PFResidual
00504                 setattr(process, jetCorrLabel[0]+'L1FastJet'   , ak5PFL1Fastjet.clone ( algorithm=jetCorrLabel[0]
00505                                                                                       , srcRho=cms.InputTag(newLabel('kt6'+jetCollType+'Jets'),'rho') ) )
00506                 setattr(process, jetCorrLabel[0]+'L1Offset'    , ak5PFL1Offset.clone  ( algorithm=jetCorrLabel[0] ) )
00507                 setattr(process, jetCorrLabel[0]+'L2Relative'  , ak5PFL2Relative.clone( algorithm=jetCorrLabel[0] ) )
00508                 setattr(process, jetCorrLabel[0]+'L3Absolute'  , ak5PFL3Absolute.clone( algorithm=jetCorrLabel[0] ) )
00509                 setattr(process, jetCorrLabel[0]+'L2L3Residual', ak5PFResidual.clone  ( algorithm=jetCorrLabel[0] ) )
00510                 ## combinded corrections
00511                 setattr(process, jetCorrLabel[0]+'CombinedCorrector', cms.ESProducer( 'JetCorrectionESChain'
00512                                                                                   , correctors = cms.vstring() ) )
00513 
00514                 for corrLbl in jetCorrLabel[1]:
00515                     if corrLbl != 'L1FastJet' and corrLbl != 'L1Offset' and corrLbl != 'L2Relative' and corrLbl != 'L3Absolute' and corrLbl != 'L2L3Residual':
00516                         print '========================================='
00517                         print ' Type1MET corrections are currently only  '
00518                         print ' supported for the following corrections: '
00519                         print '   - L1FastJet'
00520                         print '   - L1Offset'
00521                         print '   - L2Relative'
00522                         print '   - L3Absolute'
00523                         print '   - L2L3Residual'
00524                         print ' But given was:'
00525                         print '   -', corrLbl
00526                         print '============================================='
00527                         raise ValueError, 'unsupported JEC for TypeI MET correction: '+corrLbl
00528                     else:
00529                         getattr(process, jetCorrLabel[0]+'CombinedCorrector').correctors.append(jetCorrLabel[0]+corrLbl)
00530 
00531                 ## configuration of MET corrections
00532                 if jetCollType == 'Calo':
00533                     from JetMETCorrections.Type1MET.caloMETCorrections_cff import caloJetMETcorr,caloType1CorrectedMet,caloType1p2CorrectedMet,produceCaloMETCorrections
00534 
00535                     setattr(process,jetCorrLabel[0]+'JetMETcorr',   caloJetMETcorr.clone(srcMET       = "corMetGlobalMuons"))
00536                     setattr(process,jetCorrLabel[0]+'Type1CorMet',  caloType1CorrectedMet.clone(src   = "corMetGlobalMuons"))
00537                     setattr(process,jetCorrLabel[0]+'Type1p2CorMet',caloType1p2CorrectedMet.clone(src = "corMetGlobalMuons"))
00538 
00539                     getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).src          = cms.InputTag(jetCollection.getModuleLabel())
00540                     if ('L1FastJet' in jetCorrLabel[1] or 'L1Fastjet' in jetCorrLabel[1]):
00541                         getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).offsetCorrLabel = cms.string(jetCorrLabel[0]+'L1FastJet')
00542                     elif ('L1Offset' in jetCorrLabel[1]):
00543                         getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).offsetCorrLabel = cms.string(jetCorrLabel[0]+'L1Offset')
00544                     else:
00545                         getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).offsetCorrLabel = cms.string('')
00546                     getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).jetCorrLabel = cms.string(jetCorrLabel[0]+'CombinedCorrector')
00547 
00548                     getattr(process,jetCorrLabel[0]+'Type1CorMet'  ).srcType1Corrections = cms.VInputTag(
00549                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'type1')
00550                         )
00551 
00552                     getattr(process,jetCorrLabel[0]+'Type1p2CorMet').srcType1Corrections = cms.VInputTag(
00553                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'type1')
00554                         )
00555                     getattr(process,jetCorrLabel[0]+'Type1p2CorMet').srcUnclEnergySums = cms.VInputTag(
00556                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'type2'),
00557                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'offset'),
00558                         cms.InputTag('muonCaloMETcorr')
00559                         )
00560 
00561                     ## add MET corrections to sequence
00562                     setattr(process,'patMETs'+jetCorrLabel[0],getattr(process,'patMETs').clone(metSource = cms.InputTag(jetCorrLabel[0]+'Type1CorMet'),addMuonCorrections = False))
00563 
00564                     setattr(process,'produce'+jetCorrLabel[0]+'METCorrections',produceCaloMETCorrections.copy())
00565                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'caloJetMETcorr'),         getattr(process,jetCorrLabel[0]+'JetMETcorr'))
00566                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'caloType1CorrectedMet'),  getattr(process,jetCorrLabel[0]+'Type1CorMet'))
00567                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'caloType1p2CorrectedMet'),getattr(process,jetCorrLabel[0]+'Type1p2CorMet'))
00568 
00569                     process.patDefaultSequence.replace( getattr(process,'patMETs'+jetCorrLabel[0]),
00570                                                         getattr(process,'produce'+jetCorrLabel[0]+'METCorrections')
00571                                                         *getattr(process,'patMETs'+jetCorrLabel[0]))
00572 
00573                 elif jetCollType == 'PF':
00574                     from JetMETCorrections.Type1MET.pfMETCorrections_cff import pfCandsNotInJet,pfJetMETcorr,pfCandMETcorr,pfType1CorrectedMet,pfType1p2CorrectedMet,producePFMETCorrections
00575                     setattr(process,jetCorrLabel[0]+'CandsNotInJet',pfCandsNotInJet.clone(topCollection = jetCollection))
00576                     setattr(process,jetCorrLabel[0]+'JetMETcorr',   pfJetMETcorr.clone(src              = jetCollection))
00577                     setattr(process,jetCorrLabel[0]+'CandMETcorr',  pfCandMETcorr.clone(src             = cms.InputTag(jetCorrLabel[0]+'CandsNotInJet')))
00578                     setattr(process,jetCorrLabel[0]+'Type1CorMet',  pfType1CorrectedMet.clone())
00579                     setattr(process,jetCorrLabel[0]+'Type1p2CorMet',pfType1p2CorrectedMet.clone())
00580 
00581                     if ('L1FastJet' in jetCorrLabel[1] or 'L1Fastjet' in jetCorrLabel[1]):
00582                         getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).offsetCorrLabel = cms.string(jetCorrLabel[0]+'L1FastJet')
00583                     elif ('L1Offset' in jetCorrLabel[1]):
00584                         getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).offsetCorrLabel = cms.string(jetCorrLabel[0]+'L1Offset')
00585                     else:
00586                         getattr(process,jetCorrLabel[0]+'JetMETcorr'   ).offsetCorrLabel = cms.string('')
00587                     getattr(process,jetCorrLabel[0]+'JetMETcorr').jetCorrLabel    = cms.string(jetCorrLabel[0]+'CombinedCorrector')
00588 
00589                     getattr(process,jetCorrLabel[0]+'Type1CorMet').applyType0Corrections = cms.bool(False)
00590                     getattr(process,jetCorrLabel[0]+'Type1CorMet').srcType1Corrections = cms.VInputTag(
00591                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'type1')
00592                         )
00593                     getattr(process,jetCorrLabel[0]+'Type1p2CorMet').srcType1Corrections = cms.VInputTag(
00594                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'type1')
00595                         )
00596                     getattr(process,jetCorrLabel[0]+'Type1p2CorMet').applyType0Corrections = cms.bool(False)
00597                     getattr(process,jetCorrLabel[0]+'Type1p2CorMet').srcUnclEnergySums = cms.VInputTag(
00598                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'type2'),
00599                         cms.InputTag(jetCorrLabel[0]+'JetMETcorr', 'offset'),
00600                         cms.InputTag(jetCorrLabel[0]+'CandMETcorr')
00601                         )
00602 
00603                     ## add MET corrections to sequence
00604                     setattr(process,'patMETs'+jetCorrLabel[0],getattr(process,'patMETs').clone(metSource = cms.InputTag(jetCorrLabel[0]+'Type1CorMet'),addMuonCorrections = False))
00605 
00606                     setattr(process,'produce'+jetCorrLabel[0]+'METCorrections',producePFMETCorrections.copy())
00607                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'pfCandsNotInJet'),      getattr(process,jetCorrLabel[0]+'CandsNotInJet'))
00608                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'pfJetMETcorr'),         getattr(process,jetCorrLabel[0]+'JetMETcorr'))
00609                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'pfCandMETcorr'),        getattr(process,jetCorrLabel[0]+'CandMETcorr'))
00610                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'pfType1CorrectedMet'),  getattr(process,jetCorrLabel[0]+'Type1CorMet'))
00611                     getattr(process,'produce'+jetCorrLabel[0]+'METCorrections').replace(getattr(process,'pfType1p2CorrectedMet'),getattr(process,jetCorrLabel[0]+'Type1p2CorMet'))
00612 
00613                     process.patDefaultSequence.replace( getattr(process,'patMETs'+jetCorrLabel[0]),
00614                                                         getattr(process,'produce'+jetCorrLabel[0]+'METCorrections')
00615                                                         *getattr(process,'patMETs'+jetCorrLabel[0]))
00616 
00617         else:
00618             ## switch jetCorrFactors off
00619             l1Jets.addJetCorrFactors = False
00620 
00621 addJetCollection=AddJetCollection()
00622 


Member Data Documentation

Definition at line 238 of file jetTools.py.

tuple jetTools::AddJetCollection::_defaultParameters = dicttypes.SortedKeysDict() [static, private]

Definition at line 236 of file jetTools.py.

string jetTools::AddJetCollection::_label = 'addJetCollection' [static, private]

Definition at line 235 of file jetTools.py.

Definition at line 238 of file jetTools.py.