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  const edm::ParameterSet & p, bool apinew)
28  : G4VPhysicsConstructor(name)
29 {
30  myConfig = p;
31  if(apinew) {
32  dfactor = p.getParameter<double>("DarkMPFactor");
33  fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
34  } else {
35  // this is left for backward compatibility
36  dfactor = p.getParameter<double>("dark_factor");
37  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
38  }
39  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
40  particleDefFilePath = fp.fullPath();
42  fDecayProcess.reset(nullptr);
43  myHelper.reset(nullptr);
44 
45  edm::LogInfo("SimG4CoreCustomPhysics")
46  <<"CustomPhysicsListSS: Path for custom particle definition file: \n"
48 }
49 
51 }
52 
54  edm::LogInfo("SimG4CoreCustomPhysicsSS")
55  << "===== CustomPhysicsList::ConstructParticle ";
56  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
57 }
58 
60 
61  edm::LogInfo("SimG4CoreCustomPhysicsSS")
62  <<"CustomPhysicsListSS: adding CustomPhysics processes";
63 
64  fDecayProcess.reset(new G4Decay());
65  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
66 
67  for(auto particle : fParticleFactory.get()->GetCustomParticles()) {
68 
69  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
70  if(cp) {
71  G4ProcessManager* pmanager = particle->GetProcessManager();
72  edm::LogInfo("SimG4CoreCustomPhysics")
73  <<"CustomPhysicsListSS: " << particle->GetParticleName()
74  <<" PDGcode= " << particle->GetPDGEncoding()
75  << " Mass= " << particle->GetPDGMass()/GeV <<" GeV.";
76  if(pmanager) {
77  if(particle->GetPDGCharge() != 0.0) {
78  ph->RegisterProcess(new G4hMultipleScattering, particle);
79  ph->RegisterProcess(new G4hIonisation, particle);
80  }
81  if(fDecayProcess.get()->IsApplicable(*particle)) {
82  ph->RegisterProcess(fDecayProcess.get(), particle);
83  }
84  if(cp->GetCloud() && fHadronicInteraction &&
85  CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
86  edm::LogInfo("SimG4CoreCustomPhysics")
87  <<"CustomPhysicsListSS: " << particle->GetParticleName()
88  <<" CloudMass= " <<cp->GetCloud()->GetPDGMass()/GeV
89  <<" GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass()/GeV <<" GeV.";
90 
91  if(!myHelper.get()) { myHelper.reset(new G4ProcessHelper(myConfig, fParticleFactory.get())); }
92  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
93  }
94  if(particle->GetParticleType()=="darkpho"){
96  pmanager->AddDiscreteProcess(darkGamma);
97  }
98  }
99  }
100  }
101 }
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, bool useuni=false)