Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-07-18 00:48:13

0001 import FWCore.ParameterSet.Config as cms
0002 
0003 from Configuration.Eras.Era_Phase2C17I13M9_cff import Phase2C17I13M9
0004 process = cms.Process('testHGCalRecoLocal',Phase2C17I13M9)
0005 
0006 # import of standard configurations
0007 process.load('Configuration.StandardSequences.Services_cff')
0008 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0009 process.load('FWCore.MessageService.MessageLogger_cfi')
0010 process.load('Configuration.EventContent.EventContent_cff')
0011 process.load('SimGeneral.MixingModule.mixNoPU_cfi')
0012 process.load('Configuration.Geometry.GeometryExtended2026D110Reco_cff')
0013 process.load('Configuration.StandardSequences.MagneticField_cff')
0014 process.load('Configuration.StandardSequences.Generator_cff')
0015 process.load('IOMC.EventVertexGenerators.VtxSmearedRealistic50ns13TeVCollision_cfi')
0016 process.load('GeneratorInterface.Core.genFilterSummary_cff')
0017 process.load('Configuration.StandardSequences.SimIdeal_cff')
0018 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0019 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0020 process.load('Configuration.StandardSequences.Digi_cff')
0021 process.load('Configuration.StandardSequences.SimL1Emulator_cff')
0022 process.load('Configuration.StandardSequences.L1TrackTrigger_cff')
0023 process.load('Configuration.StandardSequences.DigiToRaw_cff')
0024 process.load('HLTrigger.Configuration.HLT_Fake2_cff')
0025 process.load('Configuration.StandardSequences.RawToDigi_cff')
0026 process.load('Configuration.StandardSequences.L1Reco_cff')
0027 process.load('Configuration.StandardSequences.Reconstruction_cff')
0028 process.load('Configuration.StandardSequences.RecoSim_cff')
0029 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0030 
0031 process.maxEvents = cms.untracked.PSet(
0032     input = cms.untracked.int32(1000)
0033 )
0034 
0035 process.MessageLogger.cerr.FwkReport.reportEvery = 5
0036 if hasattr(process,'MessageLogger'):
0037     process.MessageLogger.ValidHGCal=dict()
0038     process.MessageLogger.HGCalGeom=dict()
0039 
0040 # Input source
0041 process.source = cms.Source("EmptySource")
0042 
0043 process.options = cms.untracked.PSet(
0044     wantSummary = cms.untracked.bool(True),
0045     numberOfConcurrentRuns = cms.untracked.uint32(1),
0046     numberOfStreams = cms.untracked.uint32(0),
0047     numberOfThreads = cms.untracked.uint32(1),
0048     printDependencies = cms.untracked.bool(False),
0049     sizeOfStackForThreadsInKB = cms.optional.untracked.uint32,
0050 )
0051 
0052 # Production Info
0053 process.configurationMetadata = cms.untracked.PSet(
0054     version = cms.untracked.string(''),
0055     annotation = cms.untracked.string(''),
0056     name = cms.untracked.string('Applications')
0057 )
0058 
0059 
0060 # Output definition
0061 process.output = cms.OutputModule("PoolOutputModule",
0062     splitLevel = cms.untracked.int32(0),
0063     eventAutoFlushCompressedSize = cms.untracked.int32(5242880),
0064     outputCommands = cms.untracked.vstring(
0065         'keep *_*hbhe*_*_*',
0066     'keep *_g4SimHits_*_*',
0067 #       'keep *_mix_*_*',
0068     'keep *_*HGC*_*_*',
0069         ),
0070     fileName = cms.untracked.string('file:testHGCalSimWatcherV17.root'),
0071     dataset = cms.untracked.PSet(
0072         filterName = cms.untracked.string(''),
0073         dataTier = cms.untracked.string('GEN-SIM-DIGI-RAW-RECO')
0074     ),
0075     SelectEvents = cms.untracked.PSet(
0076         SelectEvents = cms.vstring('generation_step')
0077     )
0078 )
0079 
0080 # Additional output definition
0081 process.g4SimHits.Watchers = cms.VPSet(cms.PSet(
0082     SimG4HGCalValidation = cms.PSet(
0083         Names = cms.vstring(
0084         'HGCalEECellSensitive',  
0085                 'HGCalHESiliconCellSensitive',
0086         'HGCalHEScintillatorSensitive',
0087         ),
0088     Types          = cms.vint32(0,0,0),
0089     DetTypes       = cms.vint32(0,1,2),
0090     LabelLayerInfo = cms.string("HGCalInfoLayer"),
0091     Verbosity      = cms.untracked.int32(0),
0092     ),
0093     type = cms.string('SimG4HGCalValidation')
0094 ))
0095 
0096 # Other statements
0097 process.genstepfilter.triggerConditions=cms.vstring("generation_step")
0098 from Configuration.AlCa.GlobalTag import GlobalTag
0099 process.GlobalTag = GlobalTag(process.GlobalTag, 'auto:phase2_realistic_T21', '')
0100 
0101 process.generator = cms.EDProducer("FlatRandomPtGunProducer",
0102     PGunParameters = cms.PSet(
0103         MaxPt = cms.double(20.0),
0104         MinPt = cms.double(20.0),
0105         #PartID = cms.vint32(11), #--->electron
0106         PartID = cms.vint32(13), #--->muon
0107         #PartID = cms.vint32(211), #--->pion
0108         MaxEta = cms.double(3.0),
0109         MaxPhi = cms.double(3.14159265359),
0110         MinEta = cms.double(1.2),
0111         MinPhi = cms.double(-3.14159265359)
0112     ),
0113     Verbosity = cms.untracked.int32(0),
0114     psethack = cms.string('single muon pt 35'),
0115     AddAntiParticle = cms.bool(False),
0116     firstRun = cms.untracked.uint32(1)
0117 )
0118 
0119 
0120 #Modified to produce hgceedigis
0121 process.mix.digitizers = cms.PSet(process.theDigitizersValid)
0122 
0123 process.ProductionFilterSequence = cms.Sequence(process.generator)
0124 
0125 # Path and EndPath definitions
0126 process.generation_step = cms.Path(process.pgen)
0127 process.simulation_step = cms.Path(process.psim)
0128 process.genfiltersummary_step = cms.EndPath(process.genFilterSummary)
0129 process.digitisation_step = cms.Path(process.pdigi_valid)
0130 process.L1simulation_step = cms.Path(process.SimL1Emulator)
0131 process.L1TrackTrigger_step = cms.Path(process.L1TrackTrigger)
0132 process.digi2raw_step = cms.Path(process.DigiToRaw)
0133 process.raw2digi_step = cms.Path(process.RawToDigi)
0134 process.L1Reco_step = cms.Path(process.L1Reco)
0135 process.reconstruction_step = cms.Path(process.localreco)
0136 process.recosim_step = cms.Path(process.recosim)
0137 process.out_step = cms.EndPath(process.output)
0138 
0139 # Schedule definition
0140 process.schedule = cms.Schedule(process.generation_step,
0141                 process.simulation_step,
0142                 process.digitisation_step,
0143                                 process.L1simulation_step,
0144                                 process.L1TrackTrigger_step,
0145                                 process.digi2raw_step,
0146 #                                process.raw2digi_step,
0147 #                                process.L1Reco_step,
0148 #                                process.reconstruction_step,
0149 #                                process.recosim_step,
0150                                 process.out_step
0151                 )
0152 
0153 # filter all path with the production filter sequence
0154 for path in process.paths:
0155         if getattr(process,path)._seq is not None: getattr(process,path)._seq = process.ProductionFilterSequence * getattr(process,path)._seq