CMS 3D CMS Logo

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