File indexing completed on 2025-06-03 00:12:32
0001 import FWCore.ParameterSet.Config as cms
0002 from PhysicsTools.BPHNano.common_cff import *
0003
0004
0005
0006 KshortToPiPi = cms.EDProducer(
0007 'V0ReBuilder',
0008 V0s = cms.InputTag('slimmedKshortVertices'),
0009 trkSelection = cms.string('pt > 0.35 && abs(eta) < 3.0 && trackHighPurity()'),
0010 V0Selection = cms.string('0.3 < mass && mass < 0.7'),
0011 postVtxSelection = cms.string('0.3 < mass && mass < 0.7'
0012 '&& userFloat("sv_prob") > 0.0001'),
0013 beamSpot = cms.InputTag("offlineBeamSpot"),
0014 track_match = cms.InputTag('tracksBPH', 'SelectedTracks'),
0015 isLambda = cms.bool(False)
0016 )
0017
0018 LambdaToProtonPi = cms.EDProducer(
0019 'V0ReBuilder',
0020 V0s = cms.InputTag('slimmedLambdaVertices'),
0021 trkSelection = cms.string('pt > 0.35 && abs(eta) < 3.0 && trackHighPurity()'),
0022 V0Selection = cms.string('1 < mass && mass < 1.2'),
0023 postVtxSelection = cms.string('1 < mass && mass < 1.17'
0024 '&& userFloat("sv_prob") > 0.0001'),
0025 beamSpot = cms.InputTag("offlineBeamSpot"),
0026 track_match = cms.InputTag('tracksBPH', 'SelectedTracks'),
0027 isLambda = cms.bool(True)
0028 )
0029
0030
0031
0032
0033 KshortToPiPiTable = cms.EDProducer(
0034 'SimpleCompositeCandidateFlatTableProducer',
0035 src = cms.InputTag('KshortToPiPi','SelectedV0Collection'),
0036 cut = cms.string(""),
0037 name = cms.string("Kshort"),
0038 doc = cms.string("Kshort Variables"),
0039 singleton=cms.bool(False),
0040 extension=cms.bool(False),
0041 variables = cms.PSet(
0042
0043 CandVars,
0044
0045 chi2 = Var("userFloat('sv_chi2')", float, doc = "chi2 of fitted vertex", precision=12),
0046 svprob = Var("userFloat('sv_prob')", float, doc = "vertex probability of fitted vertex"),
0047 l_xy = Var("userFloat('l_xy')", float, doc = "post-fit vertex displacement on transverse plane"),
0048 l_xy_unc = Var("userFloat('l_xy_unc')", float, doc = "post-fit vertex uncertainty of the diplacement on the transverse plane"),
0049 prefit_mass = Var("userFloat('prefit_mass')", float, doc = "pre-fit mass of the vertex", precision=12),
0050 dca = Var("userFloat('dca')", float, doc = "DCA of B candidate wrt to beamspot", precision=12),
0051 dcaErr = Var("userFloat('dcaErr')", float, doc = "Error DCA of B candidate wrt to beamspot", precision=12),
0052 vtx_x = Var("userFloat('vtx_x')", float, doc = "x position of fitted vertex", precision=12),
0053 vtx_y = Var("userFloat('vtx_y')", float, doc = "y position of fitted vertex", precision=12),
0054 vtx_z = Var("userFloat('vtx_z')", float, doc = "z position of fitted vertex", precision=12),
0055 vtx_cxx = Var("userFloat('vtx_cxx')", float, doc = "error x of fitted vertex", precision=12),
0056 vtx_cyy = Var("userFloat('vtx_cyy')", float, doc = "error y of fitted vertex", precision=12),
0057 vtx_czz = Var("userFloat('vtx_czz')", float, doc = "error z of fitted vertex", precision=12),
0058 vtx_cyx = Var("userFloat('vtx_cyx')", float, doc = "error yx of fitted vertex", precision=12),
0059 vtx_czx = Var("userFloat('vtx_czx')", float, doc = "error zx of fitted vertex", precision=12),
0060 vtx_czy = Var("userFloat('vtx_czy')", float, doc = "error zy of fitted vertex", precision=12),
0061 fit_cos_theta_2D = Var("userFloat('fitted_cos_theta_2D')", float, doc = "cos 2D of fitted vertex wrt beamspot"),
0062
0063 fit_mass = Var("userFloat('fitted_mass')", float, doc = "post-fit mass of the vertex"),
0064 fit_massErr = Var("userFloat('massErr')", float, doc = "post-fit mass error", precision=12),
0065 fit_trk1_pt = Var("userFloat('trk1_pt')", float, doc = "post-fit pt of the leading track", precision=12),
0066 fit_trk1_eta = Var("userFloat('trk1_eta')", float, doc = "post-fit eta of the leading track", precision=12),
0067 fit_trk1_phi = Var("userFloat('trk1_phi')", float, doc = "post-fit phi of the leading track", precision=12),
0068 fit_trk2_pt = Var("userFloat('trk2_pt')", float, doc = "post-fit pt of the subleading track", precision=12),
0069 fit_trk2_eta = Var("userFloat('trk2_eta')", float, doc = "post-fit eta of the subleading track", precision=12),
0070 fit_trk2_phi = Var("userFloat('trk2_phi')", float, doc = "post-fit phi of the subleading track", precision=12),
0071 fit_pt = Var("userFloat('fitted_pt')", float, doc = "post-fit pt of the V0 candidate", precision=12),
0072 fit_eta = Var("userFloat('fitted_eta')", float, doc = "post-fit eta of the V0 candidate", precision=12),
0073 fit_phi = Var("userFloat('fitted_phi')", float, doc = "post-fit phi of the V0 candidate", precision=12),
0074
0075 trk1_idx = Var("userInt('trk1_idx')", int, doc = "leading track index to the BPH track collection"),
0076 trk2_idx = Var("userInt('trk2_idx')", int, doc = "subleading track index to the BPH track collection"),
0077 fit_trk1_charge = Var("userInt('fit_trk1_charge')", int, doc= "leading track charge"),
0078 fit_trk2_charge = Var("userInt('fit_trk2_charge')", int, doc= "subleading track charge"),
0079 )
0080 )
0081
0082 CountKshortToPiPi = cms.EDFilter("PATCandViewCountFilter",
0083 minNumber = cms.uint32(1),
0084 maxNumber = cms.uint32(999999),
0085 src = cms.InputTag('KshortToPiPi','SelectedV0Collection')
0086 )
0087
0088 LambdaToProtonPiTable = KshortToPiPiTable.clone(
0089 src = cms.InputTag('LambdaToProtonPi','SelectedV0Collection'),
0090 name = cms.string("Lambda"),
0091 doc = cms.string("Lambda Variable")
0092 )
0093
0094 CountLambdaToProtonPi = cms.EDFilter("PATCandViewCountFilter",
0095 minNumber = cms.uint32(1),
0096 maxNumber = cms.uint32(999999),
0097 src = cms.InputTag('LambdaToProtonPi','SelectedV0Collection')
0098 )
0099
0100
0101 KshortPiPiBPHMCMatch = cms.EDProducer("MCMatcher",
0102 src = KshortToPiPiTable.src,
0103 matched = cms.InputTag("finalGenParticlesBPH"),
0104 mcPdgId = cms.vint32(310),
0105 checkCharge = cms.bool(False),
0106 mcStatus = cms.vint32(2),
0107 maxDeltaR = cms.double(0.3),
0108 maxDPtRel = cms.double(1.0),
0109 resolveAmbiguities = cms.bool(True),
0110 resolveByMatchQuality = cms.bool(True),
0111 )
0112
0113 LambdaProtonPiBPHMCMatch = cms.EDProducer("MCMatcher",
0114 src = LambdaToProtonPiTable.src,
0115 matched = cms.InputTag("finalGenParticlesBPH"),
0116 mcPdgId = cms.vint32(3122),
0117 checkCharge = cms.bool(False),
0118 mcStatus = cms.vint32(2),
0119 maxDeltaR = cms.double(0.3),
0120 maxDPtRel = cms.double(1.0),
0121 resolveAmbiguities = cms.bool(True),
0122 resolveByMatchQuality = cms.bool(True),
0123 )
0124
0125
0126 KshortPiPiBPHMCTable = cms.EDProducer("CandMCMatchTableProducer",
0127 src = KshortToPiPiTable.src,
0128 mcMap = cms.InputTag("KshortPiPiBPHMCMatch"),
0129 objName = KshortToPiPiTable.name,
0130 objType = cms.string("Other"),
0131 branchName = cms.string("genPart"),
0132 docString = cms.string("MC matching to status==1 muons"),
0133 )
0134
0135
0136 LambdaProtonPiBPHMCTable = cms.EDProducer("CandMCMatchTableProducer",
0137 src = LambdaToProtonPiTable.src,
0138 mcMap = cms.InputTag("LambdaProtonPiBPHMCMatch"),
0139 objName = LambdaToProtonPiTable.name,
0140 objType = cms.string("Other"),
0141 branchName = cms.string("genPart"),
0142 docString = cms.string("MC matching to status==1 muons"),
0143 )
0144
0145 KshortToPiPiSequence = cms.Sequence( KshortToPiPi )
0146 KshortToPiPiSequenceMC = cms.Sequence( KshortToPiPi +KshortPiPiBPHMCMatch)
0147 KshortToPiPiTables = cms.Sequence( KshortToPiPiTable)
0148 KshortToPiPiTablesMC = cms.Sequence( KshortToPiPiTable+KshortPiPiBPHMCTable)
0149
0150
0151 LambdaToProtonPiSequence = cms.Sequence( LambdaToProtonPi )
0152 LambdaToProtonPiSequenceMC = cms.Sequence( LambdaToProtonPi + LambdaProtonPiBPHMCMatch)
0153 LambdaToProtonPiTables = cms.Sequence( LambdaToProtonPiTable)
0154 LambdaToProtonPiTablesMC = cms.Sequence( LambdaToProtonPiTable+LambdaProtonPiBPHMCTable)
0155
0156