Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-06-07 22:23:02

0001 ###############################################################################
0002 # Way to use this:
0003 #   cmsRun protoHGCalSimWatcher_cfg.py geometry=D77
0004 #
0005 #   Options for geometry D49, D68, D77, D83, D84, D88, D92
0006 #
0007 ###############################################################################
0008 import FWCore.ParameterSet.Config as cms
0009 import os, sys, imp, re
0010 import FWCore.ParameterSet.VarParsing as VarParsing
0011 
0012 ####################################################################
0013 ### SETUP OPTIONS
0014 options = VarParsing.VarParsing('standard')
0015 options.register('geometry',
0016                  "D88",
0017                   VarParsing.VarParsing.multiplicity.singleton,
0018                   VarParsing.VarParsing.varType.string,
0019                   "geometry of operations: D49, D68, D84, D77, D83, D88, D92")
0020 
0021 ### get and parse the command line arguments
0022 options.parseArguments()
0023 
0024 print(options)
0025 
0026 ####################################################################
0027 # Use the options
0028 
0029 if (options.geometry == "D49"):
0030     from Configuration.Eras.Era_Phase2C9_cff import Phase2C9
0031     process = cms.Process('PROD',Phase2C9)
0032     process.load('Configuration.Geometry.GeometryExtended2026D49_cff')
0033     process.load('Configuration.Geometry.GeometryExtended2026D49Reco_cff')
0034     fileCheck = 'testHGCalSimWatcherV11.root'
0035     runMode = 0
0036 elif (options.geometry == "D68"):
0037     from Configuration.Eras.Era_Phase2C12_cff import Phase2C12
0038     process = cms.Process('PROD',Phase2C12)
0039     process.load('Configuration.Geometry.GeometryExtended2026D68_cff')
0040     process.load('Configuration.Geometry.GeometryExtended2026D68Reco_cff')
0041     fileCheck = 'testHGCalSimWatcherV12.root'
0042     runMode = 0
0043 elif (options.geometry == "D83"):
0044     from Configuration.Eras.Era_Phase2C11M9_cff import Phase2C11M9
0045     process = cms.Process('PROD',Phase2C11M9)
0046     process.load('Configuration.Geometry.GeometryExtended2026D83_cff')
0047     process.load('Configuration.Geometry.GeometryExtended2026D83Reco_cff')
0048     fileCheck = 'testHGCalSimWatcherV15.root'
0049     runMode = 1
0050 elif (options.geometry == "D84"):
0051     from Configuration.Eras.Era_Phase2C11_cff import Phase2C11
0052     process = cms.Process('PROD',Phase2C11)
0053     process.load('Configuration.Geometry.GeometryExtended2026D84_cff')
0054     process.load('Configuration.Geometry.GeometryExtended2026D84Reco_cff')
0055     fileCheck = 'testHGCalSimWatcherV13.root'
0056     runMode = 0
0057 elif (options.geometry == "D88"):
0058     from Configuration.Eras.Era_Phase2C11_cff import Phase2C11
0059     process = cms.Process('PROD',Phase2C11)
0060     process.load('Configuration.Geometry.GeometryExtended2026D88_cff')
0061     process.load('Configuration.Geometry.GeometryExtended2026D88Reco_cff')
0062     fileCheck = 'testHGCalSimWatcherV16.root'
0063     runMode = 1
0064 elif (options.geometry == "D92"):
0065     from Configuration.Eras.Era_Phase2C11_cff import Phase2C11
0066     process = cms.Process('PROD',Phase2C11)
0067     process.load('Configuration.Geometry.GeometryExtended2026D92_cff')
0068     process.load('Configuration.Geometry.GeometryExtended2026D92Reco_cff')
0069     fileCheck = 'testHGCalSimWatcherV17.root'
0070     runMode = 1
0071 else:
0072     from Configuration.Eras.Era_Phase2C11M9_cff import Phase2C11M9
0073     process = cms.Process('PROD',Phase2C11M9)
0074     process.load('Configuration.Geometry.GeometryExtended2026D77_cff')
0075     process.load('Configuration.Geometry.GeometryExtended2026D77Reco_cff')
0076     fileCheck = 'testHGCalSimWatcherV14.root'
0077     runMode = 0
0078 
0079 # import of standard configurations
0080 process.load('Configuration.StandardSequences.Services_cff')
0081 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0082 process.load('FWCore.MessageService.MessageLogger_cfi')
0083 process.load('Configuration.EventContent.EventContent_cff')
0084 process.load('SimGeneral.MixingModule.mixNoPU_cfi')
0085 process.load('Configuration.StandardSequences.MagneticField_cff')
0086 process.load('Configuration.StandardSequences.Generator_cff')
0087 process.load('IOMC.EventVertexGenerators.VtxSmearedRealistic50ns13TeVCollision_cfi')
0088 process.load('GeneratorInterface.Core.genFilterSummary_cff')
0089 process.load('Configuration.StandardSequences.SimIdeal_cff')
0090 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0091 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0092 process.load('Configuration.StandardSequences.Digi_cff')
0093 process.load('Configuration.StandardSequences.SimL1Emulator_cff')
0094 process.load('Configuration.StandardSequences.L1TrackTrigger_cff')
0095 process.load('Configuration.StandardSequences.DigiToRaw_cff')
0096 process.load('HLTrigger.Configuration.HLT_Fake2_cff')
0097 process.load('Configuration.StandardSequences.RawToDigi_cff')
0098 process.load('Configuration.StandardSequences.L1Reco_cff')
0099 process.load('Configuration.StandardSequences.Reconstruction_cff')
0100 process.load('Configuration.StandardSequences.RecoSim_cff')
0101 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0102 
0103 process.maxEvents = cms.untracked.PSet(
0104     input = cms.untracked.int32(1000)
0105 )
0106 
0107 process.MessageLogger.cerr.FwkReport.reportEvery = 5
0108 if hasattr(process,'MessageLogger'):
0109     process.MessageLogger.ValidHGCal=dict()
0110     process.MessageLogger.HGCalGeom=dict()
0111 
0112 # Input source
0113 process.source = cms.Source("EmptySource")
0114 
0115 process.options = cms.untracked.PSet(
0116     wantSummary = cms.untracked.bool(True),
0117     numberOfConcurrentRuns = cms.untracked.uint32(1),
0118     numberOfStreams = cms.untracked.uint32(0),
0119     numberOfThreads = cms.untracked.uint32(1),
0120     printDependencies = cms.untracked.bool(False),
0121     sizeOfStackForThreadsInKB = cms.optional.untracked.uint32,
0122 )
0123 
0124 # Production Info
0125 process.configurationMetadata = cms.untracked.PSet(
0126     version = cms.untracked.string(''),
0127     annotation = cms.untracked.string(''),
0128     name = cms.untracked.string('Applications')
0129 )
0130 
0131 # Output definition
0132 process.output = cms.OutputModule("PoolOutputModule",
0133     splitLevel = cms.untracked.int32(0),
0134     eventAutoFlushCompressedSize = cms.untracked.int32(5242880),
0135     outputCommands = cms.untracked.vstring(
0136         'keep *_*hbhe*_*_*',
0137     'keep *_g4SimHits_*_*',
0138     'keep *_*HGC*_*_*',
0139         ),
0140     fileName = cms.untracked.string(fileCheck),
0141     dataset = cms.untracked.PSet(
0142         filterName = cms.untracked.string(''),
0143         dataTier = cms.untracked.string('GEN-SIM-DIGI-RAW-RECO')
0144     ),
0145     SelectEvents = cms.untracked.PSet(
0146         SelectEvents = cms.vstring('generation_step')
0147     )
0148 )
0149 
0150 # Additional output definition
0151 if (runMode == 1):
0152     process.g4SimHits.Watchers = cms.VPSet(cms.PSet(
0153         SimG4HGCalValidation = cms.PSet(
0154             Names = cms.vstring(
0155                 'HGCalEECellSensitive',  
0156                 'HGCalHESiliconCellSensitive',
0157                 'HGCalHEScintillatorSensitive',
0158             ),
0159             Types          = cms.vint32(0,0,0),
0160             DetTypes       = cms.vint32(0,1,2),
0161             LabelLayerInfo = cms.string("HGCalInfoLayer"),
0162             Verbosity      = cms.untracked.int32(0),
0163         ),
0164         type = cms.string('SimG4HGCalValidation')
0165     ))
0166 else:
0167     process.g4SimHits.Watchers = cms.VPSet(cms.PSet(
0168         SimG4HGCalValidation = cms.PSet(
0169             Names = cms.vstring(
0170                 'HGCalEESensitive',  
0171                 'HGCalHESensitive',
0172                 'HGCalHESiliconSensitive',
0173                 'HGCalHEScintillatorSensitive',
0174             ),
0175             Types          = cms.vint32(0,0,0,0),
0176             DetTypes       = cms.vint32(0,1,1,2),
0177             LabelLayerInfo = cms.string("HGCalInfoLayer"),
0178             Verbosity      = cms.untracked.int32(0),
0179         ),
0180         type = cms.string('SimG4HGCalValidation')
0181     ))
0182 
0183 # Other statements
0184 process.genstepfilter.triggerConditions=cms.vstring("generation_step")
0185 from Configuration.AlCa.GlobalTag import GlobalTag
0186 process.GlobalTag = GlobalTag(process.GlobalTag, 'auto:phase2_realistic', '')
0187 
0188 process.generator = cms.EDProducer("FlatRandomPtGunProducer",
0189     PGunParameters = cms.PSet(
0190         MaxPt = cms.double(20.0),
0191         MinPt = cms.double(20.0),
0192         PartID = cms.vint32(13), #--->muon
0193         MaxEta = cms.double(3.0),
0194         MaxPhi = cms.double(3.14159265359),
0195         MinEta = cms.double(1.2),
0196         MinPhi = cms.double(-3.14159265359)
0197     ),
0198     Verbosity = cms.untracked.int32(0),
0199     psethack = cms.string('single muon pt 20'),
0200     AddAntiParticle = cms.bool(False),
0201     firstRun = cms.untracked.uint32(1)
0202 )
0203 
0204 
0205 #Modified to produce hgceedigis
0206 process.mix.digitizers = cms.PSet(process.theDigitizersValid)
0207 
0208 process.ProductionFilterSequence = cms.Sequence(process.generator)
0209 
0210 # Path and EndPath definitions
0211 process.generation_step = cms.Path(process.pgen)
0212 process.simulation_step = cms.Path(process.psim)
0213 process.genfiltersummary_step = cms.EndPath(process.genFilterSummary)
0214 process.digitisation_step = cms.Path(process.pdigi_valid)
0215 process.L1simulation_step = cms.Path(process.SimL1Emulator)
0216 process.L1TrackTrigger_step = cms.Path(process.L1TrackTrigger)
0217 process.digi2raw_step = cms.Path(process.DigiToRaw)
0218 process.raw2digi_step = cms.Path(process.RawToDigi)
0219 process.L1Reco_step = cms.Path(process.L1Reco)
0220 process.reconstruction_step = cms.Path(process.localreco)
0221 process.recosim_step = cms.Path(process.recosim)
0222 process.out_step = cms.EndPath(process.output)
0223 
0224 # Schedule definition
0225 process.schedule = cms.Schedule(process.generation_step,
0226                 process.simulation_step,
0227                 process.digitisation_step,
0228                 process.L1simulation_step,
0229                                 process.L1TrackTrigger_step,
0230                 process.digi2raw_step,
0231 #               process.raw2digi_step,
0232 #               process.L1Reco_step,
0233 #               process.reconstruction_step,
0234 #                               process.recosim_step,
0235                 process.out_step
0236                 )
0237 
0238 # filter all path with the production filter sequence
0239 for path in process.paths:
0240         if getattr(process,path)._seq is not None: getattr(process,path)._seq = process.ProductionFilterSequence * getattr(process,path)._seq