Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:19:27

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 #--------------------------------------------------------------------------------
0004 # select collection of "good" collision vertices
0005 
0006 selectedVerticesForPFMEtCorrType0 = cms.EDFilter("VertexSelector",
0007     src = cms.InputTag('offlinePrimaryVertices'),
0008     cut = cms.string("isValid & ndof >= 4 & chi2 > 0 & tracksSize > 0 & abs(z) < 24 & abs(position.Rho) < 2."),
0009     filter = cms.bool(False)                                          
0010 )
0011 
0012 selectedPrimaryVertexHighestPtTrackSumForPFMEtCorrType0 = cms.EDFilter("PATSingleVertexSelector",
0013     mode = cms.string('firstVertex'),
0014     vertices = cms.InputTag('selectedVerticesForPFMEtCorrType0'),
0015     filter = cms.bool(False)                                                    
0016 )
0017 #--------------------------------------------------------------------------------
0018 
0019 #--------------------------------------------------------------------------------
0020 # association of PFCandidates to vertices
0021 
0022 from RecoParticleFlow.PFTracking.particleFlowDisplacedVertex_cfi import particleFlowDisplacedVertex
0023 from TrackingTools.TransientTrack.TransientTrackBuilder_cfi import *
0024 
0025 from CommonTools.RecoUtils.pfcand_assomap_cfi import PFCandAssoMap as _PFCandAssoMap
0026 pfCandidateToVertexAssociation = _PFCandAssoMap.clone(
0027     PFCandidateCollection = cms.InputTag('particleFlow'),
0028     UseBeamSpotCompatibility = cms.untracked.bool(True),
0029     ignoreMissingCollection = cms.bool(True)
0030 )
0031 #--------------------------------------------------------------------------------
0032 
0033 #--------------------------------------------------------------------------------
0034 # produce Type 0 MET corrections
0035 
0036 pfMETcorrType0 = cms.EDProducer("Type0PFMETcorrInputProducer",
0037     srcPFCandidateToVertexAssociations = cms.InputTag('pfCandidateToVertexAssociation'),
0038     srcHardScatterVertex = cms.InputTag('selectedPrimaryVertexHighestPtTrackSumForPFMEtCorrType0'),
0039     correction = cms.PSet(
0040       # RunI correction
0041       #  formula = cms.string("-([0] + [1]*x)*(1.0 + TMath::Erf(-[2]*TMath::Power(x, [3])))"),
0042       #  par0 = cms.double(0.),
0043       #  par1 = cms.double(-0.703151),
0044       #  par2 = cms.double(0.0303531),
0045       #  par3 = cms.double(0.909209)  
0046         formula = cms.string("(x<35)?(-( [0]+x*[1]+pow(x, 2)*[2]+pow(x, 3)*[3] )):(-( [0]+35*[1]+pow(35, 2)*[2]+pow(35, 3)*[3] ))"),
0047         par0 = cms.double(-1.81414e-01),
0048         par1 = cms.double(-4.76934e-01),
0049         par2 = cms.double(8.63564e-03),
0050         par3 = cms.double(-4.94181e-05)      
0051     ),
0052     minDz = cms.double(0.2) # [cm], minimum distance required between pile-up vertices and "hard scatter" vertex
0053 )   
0054 #--------------------------------------------------------------------------------
0055 
0056 type0PFMEtCorrectionPFCandToVertexAssociationTask = cms.Task(
0057     selectedVerticesForPFMEtCorrType0,
0058     selectedPrimaryVertexHighestPtTrackSumForPFMEtCorrType0,
0059     particleFlowDisplacedVertex,
0060     pfCandidateToVertexAssociation
0061 )
0062 
0063 type0PFMEtCorrectionPFCandToVertexAssociation = cms.Sequence(
0064     type0PFMEtCorrectionPFCandToVertexAssociationTask
0065 )
0066 
0067 type0PFMEtCorrectionPFCandToVertexAssociationForValidation = cms.Sequence(
0068     type0PFMEtCorrectionPFCandToVertexAssociationTask
0069 )
0070 
0071 type0PFMEtCorrectionPFCandToVertexAssociationForValidationMiniAOD = cms.Sequence(
0072     type0PFMEtCorrectionPFCandToVertexAssociationTask
0073 )
0074 
0075 type0PFMEtCorrectionTask = cms.Task(
0076     type0PFMEtCorrectionPFCandToVertexAssociationTask,
0077     pfMETcorrType0
0078 )
0079 
0080 type0PFMEtCorrection = cms.Sequence(
0081     type0PFMEtCorrectionTask
0082 )