src
SimG4Core
CustomPhysics
src
CustomPhysics.cc
Go to the documentation of this file.
1
#include "
SimG4Core/CustomPhysics/interface/CustomPhysics.h
"
2
#include "
SimG4Core/CustomPhysics/interface/CustomPhysicsList.h
"
3
#include "
SimG4Core/CustomPhysics/interface/CustomPhysicsListSS.h
"
4
#include "
SimG4Core/PhysicsLists/interface/CMSEmStandardPhysics.h
"
5
#include "
SimG4Core/PhysicsLists/interface/CMSHadronPhysicsFTFP_BERT.h
"
6
#include "
FWCore/MessageLogger/interface/MessageLogger.h
"
7
#include "
SimG4Core/CustomPhysics/interface/APrimePhysics.h
"
8
9
#include "G4DecayPhysics.hh"
10
#include "G4EmExtraPhysics.hh"
11
#include "G4IonPhysics.hh"
12
#include "G4StoppingPhysics.hh"
13
#include "G4HadronElasticPhysics.hh"
14
#include "G4NeutronTrackingCut.hh"
15
16
#include <CLHEP/Units/SystemOfUnits.h>
17
using
CLHEP::ns;
18
19
CustomPhysics::CustomPhysics
(
const
edm::ParameterSet
&
p
) :
PhysicsList
(
p
) {
20
int
ver =
p
.getUntrackedParameter<
int
>(
"Verbosity"
, 0);
21
bool
tracking
=
p
.getParameter<
bool
>(
"TrackingCut"
);
22
bool
ssPhys =
p
.getUntrackedParameter<
bool
>(
"ExoticaPhysicsSS"
,
false
);
23
bool
dbrem =
p
.getUntrackedParameter<
bool
>(
"DBrem"
,
false
);
24
double
timeLimit =
p
.getParameter<
double
>(
"MaxTrackTime"
) * ns;
25
edm::LogInfo
(
"PhysicsList"
) <<
"You are using the simulation engine: "
26
<<
"FTFP_BERT_EMM for regular particles \n"
27
<<
"CustomPhysicsList "
<< ssPhys <<
" for exotics; "
28
<<
" tracking cut "
<<
tracking
<<
" t(ns)= "
<< timeLimit / ns;
29
// EM Physics
30
RegisterPhysics(
new
CMSEmStandardPhysics
(ver,
p
));
31
32
// Synchroton Radiation & GN Physics
33
RegisterPhysics(
new
G4EmExtraPhysics(ver));
34
35
// Decays
36
RegisterPhysics(
new
G4DecayPhysics(ver));
37
38
// Hadron Elastic scattering
39
RegisterPhysics(
new
G4HadronElasticPhysics(ver));
40
41
// Hadron Physics
42
RegisterPhysics(
new
CMSHadronPhysicsFTFP_BERT
(ver));
43
44
// Stopping Physics
45
RegisterPhysics(
new
G4StoppingPhysics(ver));
46
47
// Ion Physics
48
RegisterPhysics(
new
G4IonPhysics(ver));
49
50
// Neutron tracking cut
51
if
(
tracking
) {
52
G4NeutronTrackingCut* ncut =
new
G4NeutronTrackingCut(ver);
53
ncut->SetTimeLimit(timeLimit);
54
RegisterPhysics(ncut);
55
}
56
57
// Custom Physics
58
if
(dbrem) {
59
RegisterPhysics(
new
APrimePhysics
(
p
.getUntrackedParameter<
double
>(
"DBremMass"
),
60
p
.getUntrackedParameter<
std::string
>(
"DBremScaleFile"
),
61
p
.getUntrackedParameter<
double
>(
"DBremBiasFactor"
)));
62
}
else
if
(ssPhys) {
63
RegisterPhysics(
new
CustomPhysicsListSS
(
"custom"
,
p
));
64
}
else
{
65
RegisterPhysics(
new
CustomPhysicsList
(
"custom"
,
p
));
66
}
67
}
CMSHadronPhysicsFTFP_BERT
Definition:
CMSHadronPhysicsFTFP_BERT.h:23
MessageLogger.h
CustomPhysicsListSS
Definition:
CustomPhysicsListSS.h:11
CustomPhysicsList.h
APrimePhysics
Definition:
APrimePhysics.h:7
CMSEmStandardPhysics
Definition:
CMSEmStandardPhysics.h:17
APrimePhysics.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition:
AlCaHLTBitMon_QueryRunRegistry.py:256
PhysicsList
Definition:
PhysicsList.h:7
CustomPhysicsList
Definition:
CustomPhysicsList.h:12
CustomPhysics::CustomPhysics
CustomPhysics(const edm::ParameterSet &p)
Definition:
CustomPhysics.cc:19
CustomPhysicsListSS.h
CMSEmStandardPhysics.h
edm::LogInfo
Log< level::Info, false > LogInfo
Definition:
MessageLogger.h:131
edm::ParameterSet
Definition:
ParameterSet.h:48
CMSHadronPhysicsFTFP_BERT.h
AlCaHLTBitMon_ParallelJobs.p
def p
Definition:
AlCaHLTBitMon_ParallelJobs.py:153
CustomPhysics.h
tracking
Definition:
TempMeasurements.h:8
Generated for CMSSW Reference Manual by
1.8.14