1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
|
import FWCore.ParameterSet.Config as cms
from Configuration.Eras.Modifier_h2tb_cff import h2tb
process = cms.Process("PROD", h2tb)
process.load('FWCore.MessageService.MessageLogger_cfi')
process.load("SimGeneral.HepPDTESSource.pdt_cfi")
process.load("Geometry.EcalCommonData.ecalSimulationParameters_cff")
process.load("Geometry.HcalTestBeamData.hcalDDDSimConstants_cff")
process.load("Geometry.HcalTestBeamData.TB2007TestGeometryXML_cfi")
process.load('GeneratorInterface.Core.generatorSmeared_cfi')
process.load("Configuration.EventContent.EventContent_cff")
process.load("SimG4Core.Application.g4SimHits_cfi")
process.MessageLogger.cerr.FwkReport.reportEvery = 5
if hasattr(process,'MessageLogger'):
process.MessageLogger.HCalGeom=dict()
process.MessageLogger.EcalGeom=dict()
process.MessageLogger.VertexGenerator=dict()
process.Timing = cms.Service("Timing")
process.load("IOMC.RandomEngine.IOMC_cff")
process.RandomNumberGeneratorService.generator.initialSeed = 456789
process.RandomNumberGeneratorService.g4SimHits.initialSeed = 9876
process.RandomNumberGeneratorService.VtxSmeared.initialSeed = 123456789
beamPosition = -800.0
process.common_beam_direction_parameters = cms.PSet(
MaxEta = cms.double(0.5655),
MinEta = cms.double(0.5655),
MaxPhi = cms.double(-0.1309),
MinPhi = cms.double(-0.1309),
BeamPosition = cms.double(beamPosition)
)
process.source = cms.Source("EmptySource")
process.generator = cms.EDProducer("FlatRandomEGunProducer",
PGunParameters = cms.PSet(
process.common_beam_direction_parameters,
MinE = cms.double(9.99),
MaxE = cms.double(10.01),
PartID = cms.vint32(211)
),
Verbosity = cms.untracked.int32(0),
AddAntiParticle = cms.bool(False),
firstRun = cms.untracked.uint32(1)
)
from IOMC.EventVertexGenerators.VtxSmearedParameters_cfi import *
process.VtxSmeared = cms.EDProducer("BeamProfileVtxGenerator",
process.common_beam_direction_parameters,
VtxSmearedCommon,
BeamMeanX = cms.double(0.0),
BeamMeanY = cms.double(0.0),
BeamSigmaX = cms.double(0.0001),
BeamSigmaY = cms.double(0.0001),
Psi = cms.double(999.9),
GaussianProfile = cms.bool(False),
BinX = cms.int32(50),
BinY = cms.int32(50),
File = cms.string('beam.profile'),
UseFile = cms.bool(False),
TimeOffset = cms.double(0.)
)
process.maxEvents = cms.untracked.PSet(
input = cms.untracked.int32(1)
)
process.p1 = cms.Path(process.generator*process.VtxSmeared*process.generatorSmeared*process.g4SimHits)
process.g4SimHits.NonBeamEvent = True
process.g4SimHits.UseMagneticField = False
process.g4SimHits.Physics.type = 'SimG4Core/Physics/QGSP_FTFP_BERT_EML'
process.g4SimHits.Physics.Region = 'HcalRegion'
process.g4SimHits.Physics.DefaultCutValue = 1.
process.g4SimHits.StackingAction.KillGamma = False
process.g4SimHits.CaloSD.BeamPosition = beamPosition
process.g4SimHits.ECalSD.UseBirkLaw = True
process.g4SimHits.ECalSD.BirkL3Parametrization = True
process.g4SimHits.ECalSD.BirkC1 = 0.033
process.g4SimHits.ECalSD.BirkC2 = 0.0
process.g4SimHits.ECalSD.SlopeLightYield = 0.05
process.g4SimHits.HCalSD.UseBirkLaw = True
process.g4SimHits.HCalSD.BirkC1 = 0.0052
process.g4SimHits.HCalSD.BirkC2 = 0.142
process.g4SimHits.HCalSD.BirkC3 = 1.75
process.g4SimHits.HCalSD.UseLayerWt = False
process.g4SimHits.HCalSD.WtFile = ' '
process.g4SimHits.HCalSD.UseShowerLibrary = False
process.g4SimHits.HCalSD.TestNumberingScheme = True
process.g4SimHits.HCalSD.UseHF = False
process.g4SimHits.HCalSD.ForTBHCAL = True
process.g4SimHits.HCalSD.ForTBH2 = True
process.g4SimHits.OnlySDs = ['CaloTrkProcessing',
'HcalTB06BeamDetector',
'EcalSensitiveDetector',
'HcalSensitiveDetector']
process.g4SimHits.CaloTrkProcessing.TestBeam = True
# Geant4 geometry check
process.g4SimHits.G4CheckOverlap.OutputBaseName = cms.string("hctb07")
process.g4SimHits.G4CheckOverlap.OverlapFlag = cms.bool(True)
process.g4SimHits.G4CheckOverlap.Tolerance = cms.double(0.1)
process.g4SimHits.G4CheckOverlap.Resolution = cms.int32(10000)
process.g4SimHits.G4CheckOverlap.Depth = cms.int32(-1)
process.g4SimHits.G4CheckOverlap.RegionFlag = cms.bool(False)
process.g4SimHits.G4CheckOverlap.NodeNames = cms.vstring('TBHCal')
process.g4SimHits.G4CheckOverlap.gdmlFlag = cms.bool(False)
process.g4SimHits.G4CheckOverlap.PVname = ''
process.g4SimHits.G4CheckOverlap.LVname = ''
process.g4SimHits.FileNameField = ''
process.g4SimHits.FileNameGDML = ''
process.g4SimHits.FileNameRegions = ''
|