CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CustomPhysicsListSS.cc
Go to the documentation of this file.
1 #include <memory>
2 
9 
13 
14 #include "G4hMultipleScattering.hh"
15 #include "G4hIonisation.hh"
16 #include "G4CoulombScattering.hh"
17 #include "G4ProcessManager.hh"
18 
21 
22 using namespace CLHEP;
23 
24 G4ThreadLocal std::unique_ptr<G4ProcessHelper> CustomPhysicsListSS::myHelper;
25 
27  : G4VPhysicsConstructor(name) {
28  myConfig = p;
29  if (apinew) {
30  dfactor = p.getParameter<double>("DarkMPFactor");
31  fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
32  } else {
33  // this is left for backward compatibility
34  dfactor = p.getParameter<double>("dark_factor");
35  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
36  }
37  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
39  fParticleFactory = std::make_unique<CustomParticleFactory>();
40  myHelper.reset(nullptr);
41 
42  edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsListSS: Path for custom particle definition file: \n"
44 }
45 
47 
49  edm::LogVerbatim("SimG4CoreCustomPhysicsSS") << "===== CustomPhysicsList::ConstructParticle ";
50  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
51 }
52 
54  edm::LogVerbatim("SimG4CoreCustomPhysicsSS") << "CustomPhysicsListSS: adding CustomPhysics processes";
55 
56  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
57 
58  for (auto particle : fParticleFactory.get()->getCustomParticles()) {
59  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
60  if (cp) {
61  G4ProcessManager* pmanager = particle->GetProcessManager();
62  edm::LogVerbatim("SimG4CoreCustomPhysics")
63  << "CustomPhysicsListSS: " << particle->GetParticleName() << " PDGcode= " << particle->GetPDGEncoding()
64  << " Mass= " << particle->GetPDGMass() / GeV << " GeV.";
65  if (pmanager) {
66  if (particle->GetPDGCharge() != 0.0) {
67  ph->RegisterProcess(new G4CoulombScattering, particle);
68  ph->RegisterProcess(new G4hIonisation, particle);
69  }
70  if (cp->GetCloud() && fHadronicInteraction &&
71  (CustomPDGParser::s_isgluinoHadron(particle->GetPDGEncoding()) ||
72  (CustomPDGParser::s_isstopHadron(particle->GetPDGEncoding())) ||
73  (CustomPDGParser::s_issbottomHadron(particle->GetPDGEncoding())))) {
74  edm::LogVerbatim("SimG4CoreCustomPhysics")
75  << "CustomPhysicsListSS: " << particle->GetParticleName()
76  << " CloudMass= " << cp->GetCloud()->GetPDGMass() / GeV
77  << " GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass() / GeV << " GeV.";
78 
79  if (!myHelper.get()) {
80  myHelper = std::make_unique<G4ProcessHelper>(myConfig, fParticleFactory.get());
81  }
82  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
83  }
84  if (particle->GetParticleType() == "darkpho") {
86  pmanager->AddDiscreteProcess(darkGamma);
87  }
88  }
89  }
90  }
91 }
std::string particleDefFilePath
Log< level::Info, true > LogVerbatim
std::unique_ptr< CustomParticleFactory > fParticleFactory
void ConstructParticle() override
static bool s_isgluinoHadron(int pdg)
G4ParticleDefinition * GetCloud()
G4ParticleDefinition * GetSpectator()
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
edm::ParameterSet myConfig
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void ConstructProcess() override
static bool s_issbottomHadron(int pdg)
CustomPhysicsListSS(const std::string &name, const edm::ParameterSet &p, bool useuni=false)
std::string fullPath() const
Definition: FileInPath.cc:161
static bool s_isstopHadron(int pdg)