CMS 3D CMS Logo

CMSEmStandardPhysics.cc
Go to the documentation of this file.
3 
4 #include "G4SystemOfUnits.hh"
5 #include "G4ParticleDefinition.hh"
6 #include "G4EmParameters.hh"
7 #include "G4EmBuilder.hh"
8 
9 #include "G4ComptonScattering.hh"
10 #include "G4GammaConversion.hh"
11 #include "G4PhotoElectricEffect.hh"
12 #include "G4LivermorePhotoElectricModel.hh"
13 
14 #include "G4MscStepLimitType.hh"
15 
16 #include "G4eMultipleScattering.hh"
17 #include "G4hMultipleScattering.hh"
18 #include "G4eCoulombScatteringModel.hh"
19 #include "G4CoulombScattering.hh"
20 #include "G4WentzelVIModel.hh"
21 #include "G4UrbanMscModel.hh"
22 
23 #include "G4eIonisation.hh"
24 #include "G4eBremsstrahlung.hh"
25 #include "G4eplusAnnihilation.hh"
26 //#include "G4UAtomicDeexcitation.hh"
27 
28 #include "G4hIonisation.hh"
29 #include "G4ionIonisation.hh"
30 
31 #include "G4ParticleTable.hh"
32 #include "G4Gamma.hh"
33 #include "G4Electron.hh"
34 #include "G4Positron.hh"
35 #include "G4GenericIon.hh"
36 
37 #include "G4PhysicsListHelper.hh"
38 #include "G4BuilderType.hh"
39 #include "G4GammaGeneralProcess.hh"
40 #include "G4LossTableManager.hh"
41 
42 #include "G4RegionStore.hh"
43 #include "G4Region.hh"
44 #include <string>
45 
47  : G4VPhysicsConstructor("CMSEmStandard_emm"), verbose(ver) {
48  G4EmParameters* param = G4EmParameters::Instance();
49  param->SetDefaults();
50  param->SetVerbose(verbose);
51  param->SetApplyCuts(true);
52  param->SetStepFunction(0.8, 1 * CLHEP::mm);
53  param->SetMscRangeFactor(0.2);
54  param->SetMscStepLimitType(fMinimal);
55  SetPhysicsType(bElectromagnetic);
56  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
57  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
58  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
59  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
60  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
61  fStepLimitType = fUseSafety;
62  if (msc == "UseSafetyPlus") {
63  fStepLimitType = fUseSafetyPlus;
64  }
65  if (msc == "Minimal") {
66  fStepLimitType = fMinimal;
67  }
68 }
69 
71 
73  // minimal set of particles for EM physics
74  G4EmBuilder::ConstructMinimalEmSet();
75 }
76 
78  if (verbose > 0) {
79  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct EM Processes";
80  }
81 
82  // This EM builder takes default models of Geant4 10 EMV.
83  // Multiple scattering by Urban for all particles
84  // except e+e- below 100 MeV for which the Urban93 model is used
85  G4EmBuilder::PrepareEMPhysics();
86 
87  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
88  // processes used by several particles
89  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
90  G4NuclearStopping* pnuc(nullptr);
91 
92  // high energy limit for e+- scattering models and bremsstrahlung
93  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
94 
95  const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
96  const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
97 
98  // Add gamma EM Processes
99  G4ParticleDefinition* particle = G4Gamma::Gamma();
100 
101  G4PhotoElectricEffect* pee = new G4PhotoElectricEffect();
102  pee->SetEmModel(new G4LivermorePhotoElectricModel());
103 
104  if (G4EmParameters::Instance()->GeneralProcessActive()) {
105  G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
106  sp->AddEmProcess(pee);
107  sp->AddEmProcess(new G4ComptonScattering());
108  sp->AddEmProcess(new G4GammaConversion());
109  G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
110  ph->RegisterProcess(sp, particle);
111 
112  } else {
113  ph->RegisterProcess(pee, particle);
114  ph->RegisterProcess(new G4ComptonScattering(), particle);
115  ph->RegisterProcess(new G4GammaConversion(), particle);
116  }
117 
118  // e-
119  particle = G4Electron::Electron();
120 
121  G4eIonisation* eioni = new G4eIonisation();
122 
123  G4eMultipleScattering* msc = new G4eMultipleScattering;
124  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
125  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
126  G4UrbanMscModel* msc3 = new G4UrbanMscModel();
127  msc1->SetHighEnergyLimit(highEnergyLimit);
128  msc2->SetLowEnergyLimit(highEnergyLimit);
129  msc3->SetHighEnergyLimit(highEnergyLimit);
130  msc3->SetRangeFactor(fRangeFactor);
131  msc3->SetGeomFactor(fGeomFactor);
132  msc3->SetSafetyFactor(fSafetyFactor);
133  msc3->SetLambdaLimit(fLambdaLimit);
134  msc3->SetStepLimitType(fStepLimitType);
135  msc3->SetLocked(true);
136  msc->SetEmModel(msc1);
137  msc->SetEmModel(msc2);
138  if (nullptr != aRegion) {
139  msc->AddEmModel(-1, msc3, aRegion);
140  }
141  if (nullptr != bRegion) {
142  msc->AddEmModel(-1, msc3, bRegion);
143  }
144 
145  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
146  G4CoulombScattering* ss = new G4CoulombScattering();
147  ss->SetEmModel(ssm);
148  ss->SetMinKinEnergy(highEnergyLimit);
149  ssm->SetLowEnergyLimit(highEnergyLimit);
150  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
151 
152  ph->RegisterProcess(msc, particle);
153  ph->RegisterProcess(eioni, particle);
154  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
155  ph->RegisterProcess(ss, particle);
156 
157  // e+
158  particle = G4Positron::Positron();
159  eioni = new G4eIonisation();
160 
161  msc = new G4eMultipleScattering();
162  msc1 = new G4UrbanMscModel();
163  msc2 = new G4WentzelVIModel();
164  msc3 = new G4UrbanMscModel();
165  msc1->SetHighEnergyLimit(highEnergyLimit);
166  msc2->SetLowEnergyLimit(highEnergyLimit);
167  msc3->SetHighEnergyLimit(highEnergyLimit);
168  msc3->SetRangeFactor(fRangeFactor);
169  msc3->SetGeomFactor(fGeomFactor);
170  msc3->SetSafetyFactor(fSafetyFactor);
171  msc3->SetLambdaLimit(fLambdaLimit);
172  msc3->SetStepLimitType(fStepLimitType);
173  msc3->SetLocked(true);
174  msc->SetEmModel(msc1);
175  msc->SetEmModel(msc2);
176  if (nullptr != aRegion) {
177  msc->AddEmModel(-1, msc3, aRegion);
178  }
179  if (nullptr != bRegion) {
180  msc->AddEmModel(-1, msc3, bRegion);
181  }
182 
183  ssm = new G4eCoulombScatteringModel();
184  ss = new G4CoulombScattering();
185  ss->SetEmModel(ssm);
186  ss->SetMinKinEnergy(highEnergyLimit);
187  ssm->SetLowEnergyLimit(highEnergyLimit);
188  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
189 
190  ph->RegisterProcess(msc, particle);
191  ph->RegisterProcess(eioni, particle);
192  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
193  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
194  ph->RegisterProcess(ss, particle);
195 
196  // generic ion
197  particle = G4GenericIon::GenericIon();
198  G4ionIonisation* ionIoni = new G4ionIonisation();
199  ph->RegisterProcess(hmsc, particle);
200  ph->RegisterProcess(ionIoni, particle);
201 
202  // muons, hadrons ions
203  G4EmBuilder::ConstructCharged(hmsc, pnuc);
204 }
MessageLogger.h
CMSEmStandardPhysics::fLambdaLimit
G4double fLambdaLimit
Definition: CMSEmStandardPhysics.h:22
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CMSEmStandardPhysics::ConstructParticle
void ConstructParticle() override
Definition: CMSEmStandardPhysics.cc:72
CMSEmStandardPhysics::CMSEmStandardPhysics
CMSEmStandardPhysics(G4int ver, const edm::ParameterSet &p)
Definition: CMSEmStandardPhysics.cc:46
contentValuesCheck.ss
ss
Definition: contentValuesCheck.py:33
CMSEmStandardPhysics::~CMSEmStandardPhysics
~CMSEmStandardPhysics() override
Definition: CMSEmStandardPhysics.cc:70
CMSEmStandardPhysics.h
verbose
static constexpr int verbose
Definition: HLTExoticaSubAnalysis.cc:25
CMSEmStandardPhysics::fSafetyFactor
G4double fSafetyFactor
Definition: CMSEmStandardPhysics.h:21
CMSEmStandardPhysics::fGeomFactor
G4double fGeomFactor
Definition: CMSEmStandardPhysics.h:20
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
edm::ParameterSet
Definition: ParameterSet.h:47
CMSEmStandardPhysics::fRangeFactor
G4double fRangeFactor
Definition: CMSEmStandardPhysics.h:19
nanoDQM_cff.Electron
Electron
Definition: nanoDQM_cff.py:62
CMSEmStandardPhysics::fStepLimitType
G4MscStepLimitType fStepLimitType
Definition: CMSEmStandardPhysics.h:23
CMSEmStandardPhysics::verbose
G4int verbose
Definition: CMSEmStandardPhysics.h:24
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
CMSEmStandardPhysics::ConstructProcess
void ConstructProcess() override
Definition: CMSEmStandardPhysics.cc:77