Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-05-23 23:48:07

0001 import sys
0002 import FWCore.ParameterSet.Config as cms
0003 
0004 import FWCore.ParameterSet.VarParsing as VarParsing
0005 options = VarParsing.VarParsing()
0006 options.register("isPhase2",
0007                  False,
0008                  VarParsing.VarParsing.multiplicity.singleton,
0009                  VarParsing.VarParsing.varType.bool,
0010                  "is the test running with phase-2 geometry")
0011 options.register("maxEvents",
0012                  -1,
0013                  VarParsing.VarParsing.multiplicity.singleton,
0014                  VarParsing.VarParsing.varType.int,
0015                  "number of events to run")
0016 options.parseArguments()
0017 
0018 if (options.isPhase2):
0019      from Alignment.OfflineValidation.TkAlAllInOneTool.defaultInputFiles_cff import filesDefaultMC_TTbarPhase2RECO
0020 else:
0021      from Alignment.OfflineValidation.TkAlAllInOneTool.defaultInputFiles_cff import filesDefaultMC_TTBarPU
0022 
0023 from enum import Enum
0024 class RefitType(Enum):
0025      STANDARD = 1
0026      COMMON   = 2
0027  
0028 _isDA = True
0029 _isMC = True
0030 _allFromGT = True
0031 _applyBows = True
0032 _applyExtraConditions = True
0033 _theRefitter = RefitType.COMMON # RefitType.STANDARD (other option not involving filtering)
0034 _theTrackCollection = 'generalTracks' # FIXME: 'ALCARECOTkAlMinBias' once a sample is available
0035 
0036 ###################################################################
0037 # Set default phase-2 settings
0038 ###################################################################
0039 if(options.isPhase2):
0040      import Configuration.Geometry.defaultPhase2ConditionsEra_cff as _settings
0041      _PH2_GLOBAL_TAG, _PH2_ERA = _settings.get_era_and_conditions(_settings.DEFAULT_VERSION)
0042 
0043 ###################################################################
0044 # Set the era
0045 ###################################################################
0046 if(options.isPhase2):
0047      process = cms.Process("Demo", _PH2_ERA)
0048 else:
0049      from Configuration.Eras.Era_Run3_cff import Run3
0050      process = cms.Process("Demo", Run3)
0051      
0052 ###################################################################
0053 # Set the process to run multi-threaded
0054 ###################################################################
0055 process.options.numberOfThreads = 8
0056 
0057 ###################################################################
0058 # Event source and run selection
0059 ###################################################################
0060 process.source = cms.Source("PoolSource",
0061                             fileNames = (filesDefaultMC_TTbarPhase2RECO if (options.isPhase2) else filesDefaultMC_TTBarPU),
0062                             duplicateCheckMode = cms.untracked.string('checkAllFilesOpened'))
0063 
0064 runboundary = 1
0065 process.source.firstRun = cms.untracked.uint32(int(runboundary))
0066 process.maxEvents = cms.untracked.PSet(input = cms.untracked.int32(options.maxEvents))  #(10 if (options.isPhase2) else 100)
0067 
0068 ###################################################################
0069 # JSON Filtering
0070 ###################################################################
0071 if _isMC:
0072      print("############ testPVValidation_cfg.py: msg%-i: This is Simulation!")
0073      runboundary = 1
0074 else:
0075      print("############ testPVValidation_cfg.py: msg%-i: This is DATA!")
0076      import FWCore.PythonUtilities.LumiList as LumiList
0077      process.source.lumisToProcess = LumiList.LumiList(filename ='None').getVLuminosityBlockRange()
0078 
0079 ###################################################################
0080 # Messages
0081 ###################################################################
0082 process.load('FWCore.MessageService.MessageLogger_cfi')   
0083 process.MessageLogger.cerr.enable = False
0084 process.MessageLogger.PrimaryVertexValidation=dict()  
0085 process.MessageLogger.SplitVertexResolution=dict()
0086 process.MessageLogger.FilterOutLowPt=dict()  
0087 process.MessageLogger.cout = cms.untracked.PSet(
0088     enable = cms.untracked.bool(True),
0089     threshold = cms.untracked.string("INFO"),
0090     default   = cms.untracked.PSet(limit = cms.untracked.int32(0)),                       
0091     FwkReport = cms.untracked.PSet(limit = cms.untracked.int32(-1),
0092                                    reportEvery = cms.untracked.int32(1 if (options.isPhase2) else 10)),                                                      
0093     PrimaryVertexValidation = cms.untracked.PSet( limit = cms.untracked.int32(-1)),
0094     SplitVertexResolution   = cms.untracked.PSet( limit = cms.untracked.int32(-1)),
0095     FilterOutLowPt          = cms.untracked.PSet( limit = cms.untracked.int32(-1)),
0096     enableStatistics = cms.untracked.bool(True)
0097     )
0098 
0099 ####################################################################
0100 # Produce the Transient Track Record in the event
0101 ####################################################################
0102 process.load("TrackingTools.TransientTrack.TransientTrackBuilder_cfi")
0103 
0104 ####################################################################
0105 # Get the Magnetic Field
0106 ####################################################################
0107 process.load('Configuration.StandardSequences.MagneticField_cff')
0108 
0109 ###################################################################
0110 # Standard loads
0111 ###################################################################
0112 if(options.isPhase2):
0113      process.load('Configuration.Geometry.GeometryExtendedRun4DefaultReco_cff')
0114 else:
0115      process.load("Configuration.Geometry.GeometryRecoDB_cff")
0116 
0117 ####################################################################
0118 # Get the BeamSpot
0119 ####################################################################
0120 process.load("RecoVertex.BeamSpotProducer.BeamSpot_cff")
0121 
0122 ####################################################################
0123 # Get the GlogalTag
0124 ####################################################################
0125 process.load("Configuration.StandardSequences.FrontierConditions_GlobalTag_cff")
0126 from Configuration.AlCa.GlobalTag import GlobalTag
0127 process.GlobalTag = GlobalTag(process.GlobalTag, (_PH2_GLOBAL_TAG if options.isPhase2 else '125X_mcRun3_2022_realistic_v3'), '')
0128 
0129 if _allFromGT:
0130      print("############ testPVValidation_cfg.py: msg%-i: All is taken from GT")
0131 else:
0132      if(options.isPhase2):
0133           print("########## overriding of phase-2 alignment conditions is not yet supported")
0134           pass
0135      else:
0136           ####################################################################
0137           # Get Alignment constants
0138           ####################################################################
0139           from CondCore.DBCommon.CondDBSetup_cfi import *
0140           process.trackerAlignment = cms.ESSource("PoolDBESSource",CondDBSetup,
0141                                                   connect = cms.string('frontier://FrontierProd/CMS_CONDITIONS'),
0142                                                   toGet = cms.VPSet(cms.PSet(record = cms.string('TrackerAlignmentRcd'),
0143                                                                              tag = cms.string('TrackerAlignment_Upgrade2017_design_v4')
0144                                                                         )
0145                                                                )
0146                                              )
0147           process.es_prefer_trackerAlignment = cms.ESPrefer("PoolDBESSource", "trackerAlignment")
0148 
0149           ####################################################################
0150           # Get APE
0151           ####################################################################
0152           process.setAPE = cms.ESSource("PoolDBESSource",CondDBSetup,
0153                                         connect = cms.string('frontier://FrontierProd/CMS_CONDITIONS'),
0154                                         toGet = cms.VPSet(cms.PSet(record = cms.string('TrackerAlignmentErrorExtendedRcd'),
0155                                                                    tag = cms.string('TrackerAlignmentErrorsExtended_Upgrade2017_design_v0')
0156                                                               )
0157                                                      )
0158                                    )
0159           process.es_prefer_setAPE = cms.ESPrefer("PoolDBESSource", "setAPE")
0160 
0161           ####################################################################
0162           # Kinks and Bows (optional)
0163           ####################################################################
0164           if _applyBows:
0165                print("############ testPVValidation_cfg.py: msg%-i: Applying TrackerSurfaceDeformations!")
0166                process.trackerBows = cms.ESSource("PoolDBESSource",CondDBSetup,
0167                                                   connect = cms.string('frontier://FrontierProd/CMS_CONDITIONS'),
0168                                                   toGet = cms.VPSet(cms.PSet(record = cms.string('TrackerSurfaceDeformationRcd'),
0169                                                                              tag = cms.string('TrackerSurfaceDeformations_zero')
0170                                                                         )
0171                                                                )
0172                                              )
0173                process.es_prefer_Bows = cms.ESPrefer("PoolDBESSource", "trackerBows")
0174           else:
0175                print("############ testPVValidation_cfg.py: msg%-i: MultiPVValidation: Not applying TrackerSurfaceDeformations!")
0176 
0177      ####################################################################
0178      # Extra corrections not included in the GT
0179      ####################################################################
0180      if _applyExtraConditions:
0181 
0182           import CalibTracker.Configuration.Common.PoolDBESSource_cfi
0183           ##### END OF EXTRA CONDITIONS
0184  
0185      else:
0186           print("############ testPVValidation_cfg.py: msg%-i: Not applying extra calibration constants!")
0187      
0188 ####################################################################
0189 # Load and Configure event selection
0190 ####################################################################
0191 process.primaryVertexFilter = cms.EDFilter("VertexSelector",
0192                                              src = cms.InputTag("offlinePrimaryVertices"),
0193                                              cut = cms.string("!isFake && ndof > 4 && abs(z) <= 24 && position.Rho <= 2"),
0194                                              filter = cms.bool(True)
0195                                              )
0196 
0197 process.noscraping = cms.EDFilter("FilterOutScraping",
0198                                   applyfilter = cms.untracked.bool(True),
0199                                   src =  cms.untracked.InputTag(_theTrackCollection),
0200                                   debugOn = cms.untracked.bool(False),
0201                                   numtrack = cms.untracked.uint32(10),
0202                                   thresh = cms.untracked.double(0.25)
0203                                   )
0204 
0205 process.noslowpt = cms.EDFilter("FilterOutLowPt",
0206                                 applyfilter = cms.untracked.bool(True),
0207                                 src =  cms.untracked.InputTag(_theTrackCollection),
0208                                 debugOn = cms.untracked.bool(False),
0209                                 numtrack = cms.untracked.uint32(0),
0210                                 thresh = cms.untracked.int32(1),
0211                                 ptmin  = cms.untracked.double(3.),
0212                                 runControl = cms.untracked.bool(True),
0213                                 runControlNumber = cms.untracked.vuint32(int(runboundary))
0214                                 )
0215 
0216 ####################################################################
0217 # BeamSpot check
0218 ####################################################################
0219 from RecoVertex.BeamSpotProducer.beamSpotCompatibilityChecker_cfi import beamSpotCompatibilityChecker
0220 process.BeamSpotChecker = beamSpotCompatibilityChecker.clone(
0221      bsFromFile = "offlineBeamSpot::RECO",  # source of the event beamspot (in the ALCARECO files)
0222      bsFromDB = "offlineBeamSpot::@currentProcess", # source of the DB beamspot (from Global Tag) NOTE: only if dbFromEvent is True!
0223      dbFromEvent = True,
0224      warningThr = 5, # significance threshold to emit a warning message
0225      errorThr = 10,  # significance threshold to abort the job
0226 )
0227 
0228 if _isMC:
0229      process.goodvertexSkim = cms.Sequence(process.noscraping)
0230 else:
0231      process.goodvertexSkim = cms.Sequence(process.primaryVertexFilter + process.noscraping + process.noslowpt)
0232 
0233 
0234 if(_theRefitter == RefitType.COMMON):
0235 
0236      print("############ testPVValidation_cfg.py: msg%-i: using the common track selection and refit sequence!")
0237      ####################################################################
0238      # Load and Configure Common Track Selection and refitting sequence
0239      ####################################################################
0240      import Alignment.CommonAlignment.tools.trackselectionRefitting as trackselRefit
0241      process.seqTrackselRefit = trackselRefit.getSequence(process, _theTrackCollection,
0242                                                           isPVValidation=True, 
0243                                                           TTRHBuilder='WithAngleAndTemplate',
0244                                                           usePixelQualityFlag=True,
0245                                                           openMassWindow=False,
0246                                                           cosmicsDecoMode=True,
0247                                                           cosmicsZeroTesla=False,
0248                                                           momentumConstraint=None,
0249                                                           cosmicTrackSplitting=False,
0250                                                           use_d0cut=False,
0251                                                           )
0252      
0253 elif (_theRefitter == RefitType.STANDARD):
0254 
0255      print("############ testPVValidation_cfg.py: msg%-i: using the standard single refit sequence!")
0256      ####################################################################
0257      # Load and Configure Measurement Tracker Event
0258      # (needed in case NavigationSchool is set != '')
0259      ####################################################################
0260      # process.load("RecoTracker.MeasurementDet.MeasurementTrackerEventProducer_cfi") 
0261      # process.MeasurementTrackerEvent.pixelClusterProducer = '_theTrackCollection'
0262      # process.MeasurementTrackerEvent.stripClusterProducer = '_theTrackCollection'
0263      # process.MeasurementTrackerEvent.inactivePixelDetectorLabels = cms.VInputTag()
0264      # process.MeasurementTrackerEvent.inactiveStripDetectorLabels = cms.VInputTag()
0265 
0266      ####################################################################
0267      # Load and Configure TrackRefitter
0268      ####################################################################
0269      process.load("RecoTracker.TrackProducer.TrackRefitters_cff")
0270      import RecoTracker.TrackProducer.TrackRefitters_cff
0271      process.FinalTrackRefitter = RecoTracker.TrackProducer.TrackRefitter_cfi.TrackRefitter.clone()
0272      process.FinalTrackRefitter.src = _theTrackCollection
0273      process.FinalTrackRefitter.TrajectoryInEvent = True
0274      process.FinalTrackRefitter.NavigationSchool = ''
0275      process.FinalTrackRefitter.TTRHBuilder = "WithAngleAndTemplate"
0276 
0277      ####################################################################
0278      # Sequence
0279      ####################################################################
0280      process.seqTrackselRefit = cms.Sequence(process.offlineBeamSpot*
0281                                              # in case NavigatioSchool is set !='' 
0282                                              #process.MeasurementTrackerEvent*
0283                                              process.FinalTrackRefitter)     
0284 
0285 ####################################################################
0286 # Output file
0287 ####################################################################
0288 process.TFileService = cms.Service("TFileService",
0289                                    fileName=cms.string("PVValidation_test_0.root")
0290                                   )
0291 
0292 ####################################################################
0293 # Imports of parameters
0294 ####################################################################
0295 from RecoVertex.PrimaryVertexProducer.OfflinePrimaryVertices_cfi import offlinePrimaryVertices
0296 ## modify the parameters which differ
0297 FilteringParams = offlinePrimaryVertices.TkFilterParameters.clone(
0298      maxNormalizedChi2 = 5.0,  # chi2ndof < 5
0299      maxD0Significance = 5.0,  # fake cut (requiring 1 PXB hit)
0300      maxEta = 5.0,             # as per recommendation in PR #18330
0301 )
0302 
0303 ## MM 04.05.2017 (use settings as in: https://github.com/cms-sw/cmssw/pull/18330)
0304 from RecoVertex.PrimaryVertexProducer.OfflinePrimaryVertices_cfi import DA_vectParameters
0305 DAClusterizationParams = DA_vectParameters.clone()
0306 
0307 GapClusterizationParams = cms.PSet(algorithm   = cms.string('gap'),
0308                                    TkGapClusParameters = cms.PSet(zSeparation = cms.double(0.2))  # 0.2 cm max separation betw. clusters
0309                                    )
0310 
0311 ####################################################################
0312 # Deterministic annealing clustering or Gap clustering
0313 ####################################################################
0314 def switchClusterizerParameters(da):
0315      if da:
0316           print("############ testPVValidation_cfg.py: msg%-i: Running DA Algorithm!")
0317           return DAClusterizationParams
0318      else:
0319           print("############ testPVValidation_cfg.py: msg%-i: Running GAP Algorithm!")
0320           return GapClusterizationParams
0321 
0322 ####################################################################
0323 # Print table of execution parameters
0324 ####################################################################
0325 from prettytable import PrettyTable
0326 x = PrettyTable()
0327 x.field_names = ["Parameter","Value"]
0328 x.add_row(["is DA",_isDA])
0329 x.add_row(["is MC",_isMC])
0330 x.add_row(["applyBows",_applyBows])
0331 x.add_row(["all from GT",_allFromGT])
0332 x.add_row(["extra conditions",_applyExtraConditions])
0333 x.add_row(["refitter type",_theRefitter])
0334 x.add_row(["track collection",_theTrackCollection])
0335 x.add_row(["GlobalTag",process.GlobalTag.globaltag.value()])
0336 x.add_row(["# events",options.maxEvents])
0337 print(x)
0338 
0339 ####################################################################
0340 # Configure the PVValidation Analyzer module
0341 ####################################################################
0342 process.PVValidation = cms.EDAnalyzer("PrimaryVertexValidation",
0343                                       TrackCollectionTag = cms.InputTag("FinalTrackRefitter"),
0344                                       VertexCollectionTag = cms.InputTag("offlinePrimaryVertices"),
0345                                       Debug = cms.bool(False),
0346                                       storeNtuple = cms.bool(False),
0347                                       useTracksFromRecoVtx = cms.bool(False),
0348                                       isLightNtuple = cms.bool(True),
0349                                       askFirstLayerHit = cms.bool(False),
0350                                       forceBeamSpot = cms.untracked.bool(False),
0351                                       probePt = cms.untracked.double(3.),
0352                                       minPt   = cms.untracked.double(1.),
0353                                       maxPt   = cms.untracked.double(30.),
0354                                       runControl = cms.untracked.bool(True),
0355                                       runControlNumber = cms.untracked.vuint32(int(runboundary)),
0356                                       TkFilterParameters = FilteringParams,
0357                                       TkClusParameters = switchClusterizerParameters(_isDA)
0358                                       )
0359 
0360 ####################################################################
0361 # Path
0362 ####################################################################
0363 process.p = cms.Path(process.goodvertexSkim*
0364                      process.seqTrackselRefit*
0365                      process.BeamSpotChecker*
0366                      process.PVValidation)
0367 
0368 ## PV refit part
0369 process.load("TrackingTools.TransientTrack.TransientTrackBuilder_cfi")
0370 process.offlinePrimaryVerticesFromRefittedTrks  = offlinePrimaryVertices.clone()
0371 process.offlinePrimaryVerticesFromRefittedTrks.TrackLabel                                       = cms.InputTag("FinalTrackRefitter")
0372 process.offlinePrimaryVerticesFromRefittedTrks.vertexCollections.maxDistanceToBeam              = 1
0373 process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxNormalizedChi2             = 20
0374 process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.minSiliconLayersWithHits      = 5
0375 process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxD0Significance             = 5.0
0376 # as it was prior to https://github.com/cms-sw/cmssw/commit/c8462ae4313b6be3bbce36e45373aa6e87253c59
0377 process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxD0Error                    = 1.0
0378 process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.maxDzError                    = 1.0
0379 process.offlinePrimaryVerticesFromRefittedTrks.TkFilterParameters.minPixelLayersWithHits        = 2
0380 
0381 ###################################################################
0382 # The trigger filter module
0383 ###################################################################
0384 from HLTrigger.HLTfilters.triggerResultsFilter_cfi import *
0385 process.HLTFilter = triggerResultsFilter.clone(
0386      triggerConditions = cms.vstring("HLT_ZeroBias_*"),
0387      #triggerConditions = cms.vstring("HLT_HT*"),
0388      hltResults = cms.InputTag( "TriggerResults", "", "HLT" ),
0389      l1tResults = cms.InputTag( "" ),
0390      throw = cms.bool(False)
0391 )
0392 
0393 ###################################################################
0394 # The analysis module
0395 ###################################################################
0396 process.trackanalysis = cms.EDAnalyzer("GeneralPurposeTrackAnalyzer",
0397                                        TkTag  = cms.InputTag('FinalTrackRefitter'),
0398                                        isCosmics = cms.bool(False))
0399 
0400 process.vertexanalysis = cms.EDAnalyzer('GeneralPurposeVertexAnalyzer',
0401                                         ndof = cms.int32(4),
0402                                         vertexLabel = cms.InputTag('offlinePrimaryVerticesFromRefittedTrks'),
0403                                         beamSpotLabel = cms.InputTag('offlineBeamSpot'),
0404                                         Xpos = cms.double(0.1),
0405                                         Ypos = cms.double(0),
0406                                         TkSizeBin = cms.int32(100),
0407                                         TkSizeMin = cms.double(499.5),
0408                                         TkSizeMax = cms.double(-0.5),
0409                                         DxyBin = cms.int32(100),
0410                                         DxyMin = cms.double(-2000),
0411                                         DxyMax = cms.double(2000),
0412                                         DzBin = cms.int32(100),
0413                                         DzMin = cms.double(-2000),
0414                                         DzMax = cms.double(2000),
0415                                         PhiBin = cms.int32(32),
0416                                         PhiBin2D = cms.int32(12),
0417                                         PhiMin = cms.double(-3.1415926535897931),
0418                                         PhiMax = cms.double(3.1415926535897931),
0419                                         EtaBin = cms.int32(26),
0420                                         EtaBin2D = cms.int32(8),
0421                                         EtaMin = cms.double(-2.7),
0422                                         EtaMax = cms.double(2.7))
0423 
0424 ###################################################################
0425 # The PV resolution module
0426 ###################################################################
0427 process.PrimaryVertexResolution = cms.EDAnalyzer('SplitVertexResolution',
0428                                                  storeNtuple         = cms.bool(True),
0429                                                  vtxCollection       = cms.InputTag("offlinePrimaryVerticesFromRefittedTrks"),
0430                                                  trackCollection     = cms.InputTag("FinalTrackRefitter"),
0431                                                  minVertexNdf        = cms.untracked.double(10.),
0432                                                  minVertexMeanWeight = cms.untracked.double(0.5),
0433                                                  runControl = cms.untracked.bool(True),
0434                                                  runControlNumber = cms.untracked.vuint32(int(runboundary))
0435                                                  )
0436 
0437 process.p2 = cms.Path(process.HLTFilter                               +
0438                       process.seqTrackselRefit                        +
0439                       process.BeamSpotChecker                         +
0440                       process.offlinePrimaryVerticesFromRefittedTrks  +
0441                       process.PrimaryVertexResolution                 +
0442                       process.trackanalysis                           +
0443                       process.vertexanalysis
0444                       )