CMS 3D CMS Logo

CustomPhysicsList.cc
Go to the documentation of this file.
1 #include <memory>
2 
9 
13 
14 #include "G4hMultipleScattering.hh"
15 #include "G4hIonisation.hh"
16 #include "G4ProcessManager.hh"
17 
22 
23 using namespace CLHEP;
24 
25 G4ThreadLocal std::unique_ptr<G4ProcessHelper> CustomPhysicsList::myHelper;
26 
28  : G4VPhysicsConstructor(name) {
29  myConfig = p;
30  if (apinew) {
31  dfactor = p.getParameter<double>("DarkMPFactor");
32  fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
33  } else {
34  // this is left for backward compatibility
35  dfactor = p.getParameter<double>("dark_factor");
36  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
37  }
38  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
39  particleDefFilePath = fp.fullPath();
40  fParticleFactory = std::make_unique<CustomParticleFactory>();
41  myHelper.reset(nullptr);
42 
43  edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsList: Path for custom particle definition file: \n"
44  << particleDefFilePath << "\n"
45  << " dark_factor= " << dfactor;
46 }
47 
49 
51  edm::LogVerbatim("SimG4CoreCustomPhysics") << "===== CustomPhysicsList::ConstructParticle ";
52  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
53 }
54 
56  edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsList: adding CustomPhysics processes "
57  << "for the list of particles";
58 
59  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
60 
61  for (auto particle : fParticleFactory.get()->getCustomParticles()) {
62  if (particle->GetParticleType() == "simp") {
63  G4ProcessManager* pmanager = particle->GetProcessManager();
64  if (pmanager) {
66  CMSQGSPSIMPBuilder* theQGSPSIMPB = new CMSQGSPSIMPBuilder();
67  theQGSPSIMPB->Build(simpInelPr);
68  pmanager->AddDiscreteProcess(simpInelPr);
69  } else
70  edm::LogVerbatim("CustomPhysics") << " No pmanager";
71  }
72 
73  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
74  if (cp) {
75  G4ProcessManager* pmanager = particle->GetProcessManager();
76  edm::LogVerbatim("SimG4CoreCustomPhysics")
77  << "CustomPhysicsList: " << particle->GetParticleName() << " PDGcode= " << particle->GetPDGEncoding()
78  << " Mass= " << particle->GetPDGMass() / GeV << " GeV.";
79  if (pmanager) {
80  if (particle->GetPDGCharge() != 0.0) {
81  ph->RegisterProcess(new G4hMultipleScattering, particle);
82  ph->RegisterProcess(new G4hIonisation, particle);
83  }
84  if (cp->GetCloud() && fHadronicInteraction &&
85  (CustomPDGParser::s_isgluinoHadron(particle->GetPDGEncoding()) ||
86  (CustomPDGParser::s_isstopHadron(particle->GetPDGEncoding())) ||
87  (CustomPDGParser::s_issbottomHadron(particle->GetPDGEncoding())))) {
88  edm::LogVerbatim("SimG4CoreCustomPhysics")
89  << "CustomPhysicsList: " << particle->GetParticleName()
90  << " CloudMass= " << cp->GetCloud()->GetPDGMass() / GeV
91  << " GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass() / GeV << " GeV.";
92 
93  if (!myHelper.get()) {
94  myHelper = std::make_unique<G4ProcessHelper>(myConfig, fParticleFactory.get());
95  }
96  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
97  }
98  if (particle->GetParticleType() == "darkpho") {
100  pmanager->AddDiscreteProcess(darkGamma);
101  }
102  }
103  }
104  }
105 }
Log< level::Info, true > LogVerbatim
CustomPhysicsList(const std::string &name, const edm::ParameterSet &p, bool useuni=false)
void Build(CMSSIMPInelasticProcess *aP)
void ConstructProcess() override
static bool s_isgluinoHadron(int pdg)
edm::ParameterSet myConfig
std::string particleDefFilePath
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
void ConstructParticle() override
std::unique_ptr< CustomParticleFactory > fParticleFactory
static bool s_issbottomHadron(int pdg)
static bool s_isstopHadron(int pdg)
~CustomPhysicsList() override