4 #include "SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
11 #include "G4hMultipleScattering.hh"
12 #include "G4hIonisation.hh"
13 #include "G4ProcessManager.hh"
15 #include "G4LeptonConstructor.hh"
16 #include "G4MesonConstructor.hh"
17 #include "G4BaryonConstructor.hh"
18 #include "G4ShortLivedConstructor.hh"
19 #include "G4IonConstructor.hh"
21 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
22 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
24 using namespace CLHEP;
28 : G4VPhysicsConstructor(name)
33 edm::LogInfo(
"CustomPhysics")<<
"Path for custom particle definition file: "
51 LogDebug(
"CustomPhysics") <<
" CustomPhysicsList: adding CustomPhysics processes";
52 aParticleIterator->reset();
54 while((*aParticleIterator)()) {
55 G4ParticleDefinition* particle = aParticleIterator->value();
58 LogDebug(
"CustomPhysics") << particle->GetParticleName()
59 <<
", "<<particle->GetPDGEncoding()
60 <<
" is Custom. Mass is "
61 <<particle->GetPDGMass()/GeV <<
" GeV.";
66 <<
" GeV. Spectator mass is "
67 <<
static_cast<CustomParticle*
>(particle)->GetSpectator()->GetPDGMass()/GeV
70 G4ProcessManager* pmanager = particle->GetProcessManager();
72 if(particle->GetPDGCharge()/eplus != 0) {
73 pmanager->AddProcess(
new G4hMultipleScattering,-1, 1, 1);
74 pmanager->AddProcess(
new G4hIonisation, -1, 2, 2);
77 if(particle->GetParticleType()==
"rhadron" ||
78 particle->GetParticleType()==
"mesonino" ||
79 particle->GetParticleType() ==
"sbaryon"){
81 pmanager->AddDiscreteProcess(
new FullModelHadronicProcess(
myHelper));
85 else LogDebug(
"CustomPhysics") <<
" No pmanager";
101 <<
" GeV. Spectator mass is "
102 <<
static_cast<CustomParticle*
>(particle)->GetSpectator()->GetPDGMass()/GeV
105 G4ProcessManager* pmanager = particle->GetProcessManager();
108 if(particle->GetPDGCharge()/eplus != 0){
109 pmanager->AddProcess(
new G4hMultipleScattering,-1, 1, 1);
110 pmanager->AddProcess(
new G4hIonisation, -1, 2, 2);
112 pmanager->AddDiscreteProcess(
new FullModelHadronicProcess(
myHelper));
114 else LogDebug(
"CustomPhysics") <<
" No pmanager";
119 G4ProcessManager* pmanager = particle->GetProcessManager();
121 if(particle->GetPDGCharge()/eplus != 0){
122 pmanager->AddProcess(
new G4hMultipleScattering,-1, 1, 1);
123 pmanager->AddProcess(
new G4hIonisation, -1, 2, 2);
125 pmanager->AddProcess(
new G4Decay, 1, -1, 3);
127 else LogDebug(
"CustomPhysics") <<
" No pmanager";
T getParameter(std::string const &) const
edm::ParameterSet myConfig
std::string particleDefFilePath
G4ParticleDefinition * GetCloud()
void setupRHadronPhycis(G4ParticleDefinition *particle)
static void loadCustomParticles(const std::string &filePath)
virtual ~CustomPhysicsList()
void setupSUSYPhycis(G4ParticleDefinition *particle)
G4ProcessHelper * myHelper
CustomPhysicsList(std::string name, const edm::ParameterSet &p)
static bool isCustomParticle(G4ParticleDefinition *particle)
std::string fullPath() const
virtual void ConstructParticle()
virtual void ConstructProcess()