CMS 3D CMS Logo

CustomPhysicsListSS.cc
Go to the documentation of this file.
7 
11 
12 #include "G4hMultipleScattering.hh"
13 #include "G4hIonisation.hh"
14 #include "G4CoulombScattering.hh"
15 #include "G4ProcessManager.hh"
16 
19 
20 using namespace CLHEP;
21 
22 G4ThreadLocal std::unique_ptr<G4ProcessHelper> CustomPhysicsListSS::myHelper;
23 
25  : G4VPhysicsConstructor(name) {
26  myConfig = p;
27  if (apinew) {
28  dfactor = p.getParameter<double>("DarkMPFactor");
29  fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
30  } else {
31  // this is left for backward compatibility
32  dfactor = p.getParameter<double>("dark_factor");
33  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
34  }
35  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
36  particleDefFilePath = fp.fullPath();
38  myHelper.reset(nullptr);
39 
40  edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsListSS: Path for custom particle definition file: \n"
42 }
43 
45 
47  edm::LogVerbatim("SimG4CoreCustomPhysicsSS") << "===== CustomPhysicsList::ConstructParticle ";
48  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
49 }
50 
52  edm::LogVerbatim("SimG4CoreCustomPhysicsSS") << "CustomPhysicsListSS: adding CustomPhysics processes";
53 
54  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
55 
56  for (auto particle : fParticleFactory.get()->GetCustomParticles()) {
57  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
58  if (cp) {
59  G4ProcessManager* pmanager = particle->GetProcessManager();
60  edm::LogVerbatim("SimG4CoreCustomPhysics")
61  << "CustomPhysicsListSS: " << particle->GetParticleName() << " PDGcode= " << particle->GetPDGEncoding()
62  << " Mass= " << particle->GetPDGMass() / GeV << " GeV.";
63  if (pmanager) {
64  if (particle->GetPDGCharge() != 0.0) {
65  ph->RegisterProcess(new G4CoulombScattering, particle);
66  ph->RegisterProcess(new G4hIonisation, particle);
67  }
68  if (cp->GetCloud() && fHadronicInteraction && CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
69  edm::LogVerbatim("SimG4CoreCustomPhysics")
70  << "CustomPhysicsListSS: " << particle->GetParticleName()
71  << " CloudMass= " << cp->GetCloud()->GetPDGMass() / GeV
72  << " GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass() / GeV << " GeV.";
73 
74  if (!myHelper.get()) {
76  }
77  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
78  }
79  if (particle->GetParticleType() == "darkpho") {
81  pmanager->AddDiscreteProcess(darkGamma);
82  }
83  }
84  }
85  }
86 }
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()
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, bool useuni=false)