CMS 3D CMS Logo

CustomPhysicsListSS.cc
Go to the documentation of this file.
7 
11 
12 #include "G4Decay.hh"
13 #include "G4hMultipleScattering.hh"
14 #include "G4hIonisation.hh"
15 #include "G4CoulombScattering.hh"
16 #include "G4ProcessManager.hh"
17 
20 
21 using namespace CLHEP;
22 
23 G4ThreadLocal std::unique_ptr<G4Decay> CustomPhysicsListSS::fDecayProcess;
24 G4ThreadLocal std::unique_ptr<G4ProcessHelper> CustomPhysicsListSS::myHelper;
25 
27  : G4VPhysicsConstructor(name)
28 {
29  myConfig = p;
30  dfactor = p.getParameter<double>("dark_factor");
31  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
32  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
33  particleDefFilePath = fp.fullPath();
34 
36  fDecayProcess.reset(nullptr);
37  myHelper.reset(nullptr);
38 
39  edm::LogInfo("SimG4CoreCustomPhysics")
40  <<"CustomPhysicsListSS: Path for custom particle definition file: \n"
42 }
43 
45 }
46 
48  edm::LogInfo("SimG4CoreCustomPhysicsSS")
49  << "===== CustomPhysicsList::ConstructParticle ";
50  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
51 }
52 
54 
55  edm::LogInfo("SimG4CoreCustomPhysicsSS")
56  <<"CustomPhysicsListSS: adding CustomPhysics processes";
57 
58  fDecayProcess.reset(new G4Decay());
59  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
60 
61  for(auto particle : fParticleFactory.get()->GetCustomParticles()) {
62 
63  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
64  if(cp) {
65  G4ProcessManager* pmanager = particle->GetProcessManager();
66  edm::LogInfo("SimG4CoreCustomPhysics")
67  <<"CustomPhysicsListSS: " << particle->GetParticleName()
68  <<" PDGcode= " << particle->GetPDGEncoding()
69  << " Mass= " << particle->GetPDGMass()/GeV <<" GeV.";
70  if(pmanager) {
71  if(particle->GetPDGCharge() != 0.0) {
72  ph->RegisterProcess(new G4hMultipleScattering, particle);
73  ph->RegisterProcess(new G4hIonisation, particle);
74  }
75  if(fDecayProcess.get()->IsApplicable(*particle)) {
76  ph->RegisterProcess(fDecayProcess.get(), particle);
77  }
78  if(cp->GetCloud() && fHadronicInteraction &&
79  CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
80  edm::LogInfo("SimG4CoreCustomPhysics")
81  <<"CustomPhysicsListSS: " << particle->GetParticleName()
82  <<" CloudMass= " <<cp->GetCloud()->GetPDGMass()/GeV
83  <<" GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass()/GeV <<" GeV.";
84 
85  if(!myHelper.get()) { myHelper.reset(new G4ProcessHelper(myConfig, fParticleFactory.get())); }
86  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
87  }
88  if(particle->GetParticleType()=="darkpho"){
90  pmanager->AddDiscreteProcess(darkGamma);
91  }
92  }
93  }
94  }
95 }
std::string particleDefFilePath
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
std::unique_ptr< CustomParticleFactory > fParticleFactory
void ConstructParticle() override
G4ParticleDefinition * GetCloud()
static G4ThreadLocal std::unique_ptr< G4Decay > fDecayProcess
G4ParticleDefinition * GetSpectator()
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
edm::ParameterSet myConfig
static bool s_isRHadron(int pdg)
void ConstructProcess() override
CustomPhysicsListSS(const std::string &name, const edm::ParameterSet &p)