Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:27:55

0001 ###############################################################################
0002 # Way to use this:
0003 #   cmsRun protoSimValid_cfg.py geometry=D88 type=hgcalBHValidation
0004 #
0005 #   Options for geometry D88, D92, D93
0006 #               type hgcalBHValidation, hgcalSiliconValidation
0007 #
0008 ###############################################################################
0009 import FWCore.ParameterSet.Config as cms
0010 import os, sys, imp, re
0011 import FWCore.ParameterSet.VarParsing as VarParsing
0012 
0013 ############################################################
0014 ### SETUP OPTIONS
0015 options = VarParsing.VarParsing('standard')
0016 options.register('geometry',
0017                  "D88",
0018                   VarParsing.VarParsing.multiplicity.singleton,
0019                   VarParsing.VarParsing.varType.string,
0020                   "geometry of operations: D88, D92, D93")
0021 options.register ('type',
0022                   "hgcalBHValidation",
0023                   VarParsing.VarParsing.multiplicity.singleton,
0024                   VarParsing.VarParsing.varType.string,
0025                   "type of operations: hgcalBHValidation, hgcalSiliconValidation")
0026 
0027 ### get and parse the command line arguments
0028 options.parseArguments()
0029 
0030 print(options)
0031 
0032 ############################################################
0033 # Use the options
0034 
0035 from Configuration.Eras.Era_Phase2C17I13M9_cff import Phase2C17I13M9
0036 process = cms.Process('PROD',Phase2C17I13M9)
0037 
0038 geomFile = "Configuration.Geometry.GeometryExtended2026" + options.geometry + "Reco_cff"
0039 if (options.type == "hgcalSiliconValidation"):
0040     fileName = "hgcSilValid" + options.geometry + ".root"
0041 else:
0042     fileName = "hgcBHValid" + options.geometry + ".root"
0043 
0044 print("Geometry file: ", geomFile)
0045 print("Output file:   ", fileName)
0046 
0047 # import of standard configurations
0048 process.load(geomFile)
0049 process.load('Configuration.StandardSequences.Services_cff')
0050 process.load('SimGeneral.HepPDTESSource.pythiapdt_cfi')
0051 process.load('FWCore.MessageService.MessageLogger_cfi')
0052 process.load('Configuration.EventContent.EventContent_cff')
0053 process.load('SimGeneral.MixingModule.mixNoPU_cfi')
0054 process.load('Configuration.StandardSequences.MagneticField_cff')
0055 process.load('Configuration.StandardSequences.Generator_cff')
0056 process.load('IOMC.EventVertexGenerators.VtxSmearedRealistic50ns13TeVCollision_cfi')
0057 process.load('GeneratorInterface.Core.genFilterSummary_cff')
0058 process.load('Configuration.StandardSequences.SimIdeal_cff')
0059 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0060 process.load('Configuration.StandardSequences.FrontierConditions_GlobalTag_cff')
0061 process.load('Configuration.StandardSequences.Digi_cff')
0062 process.load('Configuration.StandardSequences.SimL1Emulator_cff')
0063 process.load('Configuration.StandardSequences.L1TrackTrigger_cff')
0064 process.load('Configuration.StandardSequences.DigiToRaw_cff')
0065 process.load('HLTrigger.Configuration.HLT_Fake2_cff')
0066 process.load('Configuration.StandardSequences.RawToDigi_cff')
0067 process.load('Configuration.StandardSequences.L1Reco_cff')
0068 process.load('Configuration.StandardSequences.Reconstruction_cff')
0069 process.load('Configuration.StandardSequences.RecoSim_cff')
0070 process.load('Configuration.StandardSequences.EndOfProcess_cff')
0071 
0072 process.maxEvents = cms.untracked.PSet(
0073     input = cms.untracked.int32(2000)
0074 )
0075 
0076 process.MessageLogger.cerr.FwkReport.reportEvery = 5
0077 if hasattr(process,'MessageLogger'):
0078     process.MessageLogger.HGCalGeom=dict()
0079 
0080 # Input source
0081 process.source = cms.Source("EmptySource")
0082 
0083 process.generator = cms.EDProducer("FlatRandomPtGunProducer",
0084     PGunParameters = cms.PSet(
0085         MaxPt = cms.double(35.0),
0086         MinPt = cms.double(35.0),
0087         PartID = cms.vint32(13), #--->muon
0088         MinPhi = cms.double(-3.14159265359),
0089         MaxPhi = cms.double(3.14159265359),
0090         MinEta = cms.double(1.2),
0091         MaxEta = cms.double(3.0)
0092     ),
0093     Verbosity = cms.untracked.int32(0),
0094     psethack = cms.string('single muon pt 35'),
0095     AddAntiParticle = cms.bool(True),
0096     firstRun = cms.untracked.uint32(1)
0097 )
0098 
0099 process.options = cms.untracked.PSet(
0100     wantSummary = cms.untracked.bool(True)
0101 )
0102 
0103 # Production Info
0104 process.configurationMetadata = cms.untracked.PSet(
0105     version = cms.untracked.string(''),
0106     annotation = cms.untracked.string(''),
0107     name = cms.untracked.string('Applications')
0108 )
0109 
0110 # Other statements
0111 process.genstepfilter.triggerConditions=cms.vstring("generation_step")
0112 from Configuration.AlCa.GlobalTag import GlobalTag
0113 process.GlobalTag = GlobalTag(process.GlobalTag, 'auto:phase2_realistic', '')
0114 
0115 # Additional output definition
0116 process.TFileService = cms.Service("TFileService",
0117                                    fileName = cms.string(fileName),
0118                                    closeFileFast = cms.untracked.bool(True) )
0119 
0120 #Modified to produce hgceedigis
0121 process.mix.digitizers = cms.PSet(process.theDigitizersValid)
0122 process.ProductionFilterSequence = cms.Sequence(process.generator)
0123 
0124 #Following Removes Mag Field
0125 process.g4SimHits.UseMagneticField = False
0126 process.g4SimHits.Physics.bField = cms.double(0.0)
0127 
0128 # Path and EndPath definitions
0129 process.generation_step = cms.Path(process.pgen)
0130 process.simulation_step = cms.Path(process.psim)
0131 process.genfiltersummary_step = cms.EndPath(process.genFilterSummary)
0132 process.digitisation_step = cms.Path(process.pdigi_valid)
0133 process.L1simulation_step = cms.Path(process.SimL1Emulator)
0134 process.L1TrackTrigger_step = cms.Path(process.L1TrackTrigger)
0135 process.digi2raw_step = cms.Path(process.DigiToRaw)
0136 process.raw2digi_step = cms.Path(process.RawToDigi)
0137 process.L1Reco_step = cms.Path(process.L1Reco)
0138 process.reconstruction_step = cms.Path(process.localreco)
0139 process.recosim_step = cms.Path(process.recosim)
0140 
0141 if (options.type == "hgcalSiliconValidation"):
0142     process.load('Validation.HGCalValidation.hgcalSiliconValidation_cfi')
0143     process.analysis_step = cms.Path(process.hgcalSiliconAnalysisEE+process.hgcalSiliconAnalysisHEF)
0144 else:
0145     process.load('Validation.HGCalValidation.hgcalBHValidation_cfi')
0146     process.analysis_step = cms.Path(process.hgcalBHAnalysis)
0147 
0148 # Schedule definition
0149 process.schedule = cms.Schedule(process.generation_step,
0150                 process.simulation_step,
0151                 process.digitisation_step,
0152                                 process.L1simulation_step,
0153                                 process.L1TrackTrigger_step,
0154                                 process.digi2raw_step,
0155                                 process.raw2digi_step,
0156                                 process.L1Reco_step,
0157                                 process.reconstruction_step,
0158                                 process.recosim_step,
0159                                 process.analysis_step,
0160                 )
0161 
0162 # filter all path with the production filter sequence
0163 for path in process.paths:
0164         if getattr(process,path)._seq is not None: getattr(process,path)._seq = process.ProductionFilterSequence * getattr(process,path)._seq