Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-11-08 23:56:39

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 
0004 from Configuration.Eras.Era_Phase2C9_cff import Phase2C9
0005 process = cms.Process('USER',Phase2C9)
0006 
0007 # import of standard configurations
0008 process.load('Configuration.StandardSequences.Services_cff')
0009 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0010 process.load('FWCore.MessageService.MessageLogger_cfi')
0011 process.load('Configuration.EventContent.EventContent_cff')
0012 process.load('SimGeneral.MixingModule.mixNoPU_cfi')
0013 #process.load('SimGeneral.MixingModule.mix_POISSON_average_cfi')
0014 process.load('Configuration.Geometry.GeometryExtended2026D49Reco_cff')
0015 process.load('Configuration.StandardSequences.MagneticField_cff')
0016 process.load('Configuration.StandardSequences.Digi_cff')
0017 process.load('Configuration.StandardSequences.SimL1Emulator_cff')
0018 process.load('Configuration.StandardSequences.L1TrackTrigger_cff')
0019 process.load('Configuration.StandardSequences.DigiToRaw_cff')
0020 process.load('HLTrigger.Configuration.HLT_Fake2_cff')
0021 process.load('Configuration.StandardSequences.RawToDigi_cff')
0022 process.load('Configuration.StandardSequences.L1Reco_cff')
0023 process.load('Configuration.StandardSequences.Reconstruction_cff')
0024 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0025 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0026 
0027 process.maxEvents = cms.untracked.PSet(
0028     input = cms.untracked.int32(10)
0029 )
0030 
0031 process.source = cms.Source("PoolSource",
0032     fileNames = cms.untracked.vstring(
0033         '/store/relval/CMSSW_11_2_0_pre8/RelValSingleMuPt10/GEN-SIM-RECO/112X_mcRun4_realistic_v3_2026D49noPU-v1/00000/007d817e-9c59-4dec-959b-0f227942cdf0.root'
0034     )
0035 )
0036 
0037 
0038 process.options = cms.untracked.PSet(
0039 
0040 )
0041 
0042 # Production Info
0043 process.configurationMetadata = cms.untracked.PSet(
0044     annotation = cms.untracked.string('step2 nevts:10'),
0045     name = cms.untracked.string('Applications'),
0046     version = cms.untracked.string('$Revision: 1.19 $')
0047 )
0048 
0049 # MC vertice analyzer
0050 process.load("Validation.RecoVertex.mcverticesanalyzer_cfi")
0051 process.mcverticesanalyzer.pileupSummaryCollection = cms.InputTag("addPileupInfo","","HLT")
0052 
0053 # Output definition
0054 
0055 #process.FEVTDEBUGHLToutput = cms.OutputModule("PoolOutputModule",
0056 #    dataset = cms.untracked.PSet(
0057 #        dataTier = cms.untracked.string('GEN-SIM-RECO'),
0058 #        filterName = cms.untracked.string('')
0059 #    ),
0060 #    fileName = cms.untracked.string('step2_DIGI_L1_L1TrackTrigger_DIGI2RAW_HLT_RAW2DIGI_L1Reco_RECO.root'),
0061 #    outputCommands = process.FEVTDEBUGHLTEventContent.outputCommands,
0062 #    splitLevel = cms.untracked.int32(0)
0063 #)
0064 
0065 # # # -- Trajectory producer
0066 process.load("RecoTracker.TrackProducer.TrackRefitters_cff")
0067 process.TrackRefitter.src = 'generalTracks'
0068 process.TrackRefitter.NavigationSchool = ""
0069 
0070 process.ReadLocalMeasurement = cms.EDAnalyzer("Phase2PixelNtuple",
0071                                               trackProducer = cms.InputTag("generalTracks"),
0072                                               trajectoryInput = cms.InputTag('TrackRefitter::USER'),
0073                                               #verbose = cms.untracked.bool(True),
0074                                               #picky = cms.untracked.bool(False),                                             
0075                                               ### for using track hit association
0076                                               associatePixel = cms.bool(True),
0077                                               associateStrip = cms.bool(False),
0078                                               associateRecoTracks = cms.bool(False),
0079                                               ROUList = cms.vstring('TrackerHitsPixelBarrelLowTof',
0080                                                                     'TrackerHitsPixelBarrelHighTof',
0081                                                                     'TrackerHitsPixelEndcapLowTof',
0082                                                                     'TrackerHitsPixelEndcapHighTof'),
0083                                               ttrhBuilder = cms.string("WithTrackAngle"),
0084                                               usePhase2Tracker = cms.bool(True),
0085                                               pixelSimLinkSrc = cms.InputTag("simSiPixelDigis", "Pixel"),
0086                                               phase2TrackerSimLinkSrc = cms.InputTag("simSiPixelDigis", "Tracker")
0087                                           )
0088 
0089 
0090 # Additional output definition
0091 
0092 # Other statements
0093 process.mix.digitizers = cms.PSet(process.theDigitizersValid)
0094 
0095 # This pset is specific for producing simulated events for the designers of the PROC (InnerTracker)
0096 # They need pixel RecHits where the charge is stored with high-granularity and large dinamic range
0097 
0098 from Configuration.AlCa.GlobalTag import GlobalTag
0099 process.GlobalTag = GlobalTag(process.GlobalTag, 'auto:phase2_realistic_T15', '')
0100 
0101 # Path and EndPath definitions
0102 process.digitisation_step = cms.Path(process.pdigi_valid)
0103 process.L1simulation_step = cms.Path(process.SimL1Emulator)
0104 process.L1TrackTrigger_step = cms.Path(process.L1TrackTrigger)
0105 process.digi2raw_step = cms.Path(process.DigiToRaw)
0106 process.raw2digi_step = cms.Path(process.RawToDigi)
0107 process.L1Reco_step = cms.Path(process.L1Reco)
0108 process.reconstruction_step = cms.Path(process.reconstruction)
0109 process.user_step = cms.Path(process.TrackRefitter * process.ReadLocalMeasurement * process.mcverticesanalyzer)
0110 process.endjob_step = cms.EndPath(process.endOfProcess)
0111 #process.FEVTDEBUGHLToutput_step = cms.EndPath(process.FEVTDEBUGHLToutput)
0112 
0113 # Schedule definition
0114 # process.schedule imported from cff in HLTrigger.Configuration
0115 process.schedule.insert(0, process.digitisation_step)
0116 process.schedule.insert(1, process.L1simulation_step)
0117 process.schedule.insert(2, process.L1TrackTrigger_step)
0118 process.schedule.insert(3, process.digi2raw_step)
0119 process.schedule.extend([process.raw2digi_step,process.L1Reco_step,process.reconstruction_step,process.user_step,process.endjob_step])
0120 from PhysicsTools.PatAlgos.tools.helpers import associatePatAlgosToolsTask
0121 associatePatAlgosToolsTask(process)
0122 
0123 # customisation of the process.
0124 
0125 # Automatic addition of the customisation function from HLTrigger.Configuration.customizeHLTforMC
0126 from HLTrigger.Configuration.customizeHLTforMC import customizeHLTforMC 
0127 
0128 #call to customisation function customizeHLTforMC imported from HLTrigger.Configuration.customizeHLTforMC
0129 process = customizeHLTforMC(process)
0130 
0131 # End of customisation functions
0132 
0133 # Customisation from command line
0134 
0135 #Have logErrorHarvester wait for the same EDProducers to finish as those providing data for the OutputModule
0136 from FWCore.Modules.logErrorHarvester_cff import customiseLogErrorHarvesterUsingOutputCommands
0137 process = customiseLogErrorHarvesterUsingOutputCommands(process)
0138 
0139 # Add early deletion of temporary data products to reduce peak memory need
0140 from Configuration.StandardSequences.earlyDeleteSettings_cff import customiseEarlyDelete
0141 process = customiseEarlyDelete(process)
0142 # End adding early deletion
0143 process.TFileService = cms.Service('TFileService',
0144 fileName = cms.string("pixelntuple.root")
0145 )
0146