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 #include "G4HadronicProcess.hh"
16 
19 
24 
25 
26 using namespace CLHEP;
27 
28 G4ThreadLocal std::unique_ptr<G4ProcessHelper> CustomPhysicsList::myHelper;
29 
31  : G4VPhysicsConstructor(name) {
32  myConfig = p;
33  if (apinew) {
34  dfactor = p.getParameter<double>("DarkMPFactor");
35  fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
36  } else {
37  // this is left for backward compatibility
38  dfactor = p.getParameter<double>("dark_factor");
39  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
40  }
41  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
42  particleDefFilePath = fp.fullPath();
44  myHelper.reset(nullptr);
45 
46  edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsList: Path for custom particle definition file: \n"
47  << particleDefFilePath << "\n"
48  << " dark_factor= " << dfactor;
49 }
50 
52 
54  edm::LogVerbatim("SimG4CoreCustomPhysics") << "===== CustomPhysicsList::ConstructParticle ";
55  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
56 }
57 
59  edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsList: adding CustomPhysics processes "
60  << "for the list of particles";
61 
62  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
63 
64  for (auto particle : fParticleFactory.get()->GetCustomParticles()) {
65  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
66  if (cp) {
67  G4ProcessManager* pmanager = particle->GetProcessManager();
68  edm::LogVerbatim("SimG4CoreCustomPhysics")
69  << "CustomPhysicsList: " << particle->GetParticleName() << " PDGcode= " << particle->GetPDGEncoding()
70  << " Mass= " << particle->GetPDGMass() / GeV << " GeV.";
71  if (pmanager) {
72  if (particle->GetPDGCharge() != 0.0) {
73  ph->RegisterProcess(new G4hMultipleScattering, particle);
74  ph->RegisterProcess(new G4hIonisation, particle);
75  }
76  if (cp->GetCloud() && fHadronicInteraction &&
77  (CustomPDGParser::s_isgluinoHadron(particle->GetPDGEncoding()) ||
78  (CustomPDGParser::s_isstopHadron(particle->GetPDGEncoding())) ||
79  (CustomPDGParser::s_issbottomHadron(particle->GetPDGEncoding())))) {
80  edm::LogVerbatim("SimG4CoreCustomPhysics")
81  << "CustomPhysicsList: " << particle->GetParticleName()
82  << " CloudMass= " << cp->GetCloud()->GetPDGMass() / GeV
83  << " GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass() / GeV << " GeV.";
84 
85  if (!myHelper.get()) {
87  }
88  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
89  }
90  if (particle->GetParticleType() == "darkpho") {
92  pmanager->AddDiscreteProcess(darkGamma);
93  }
94  if (particle->GetParticleName() == "anti_sexaq") {
95  // here the different sexaquark interactions get defined
96  G4HadronicProcess* sqInelPr = new G4HadronicProcess();
97  CMSSQNeutronAnnih* sqModel = new CMSSQNeutronAnnih(particle->GetPDGMass() / GeV);
98  sqInelPr->RegisterMe(sqModel);
99  CMSSQInelasticCrossSection* sqInelXS = new CMSSQInelasticCrossSection(particle->GetPDGMass() / GeV);
100  sqInelPr->AddDataSet(sqInelXS);
101  pmanager->AddDiscreteProcess(sqInelPr);
102  // add also the looping needed to simulate flat interaction probability
103  CMSSQLoopProcess* sqLoopPr = new CMSSQLoopProcess();
104  pmanager->AddContinuousProcess(sqLoopPr);
105  CMSSQLoopProcessDiscr* sqLoopPrDiscr = new CMSSQLoopProcessDiscr(particle->GetPDGMass() / GeV);
106  pmanager->AddDiscreteProcess(sqLoopPrDiscr);
107  } else if (particle->GetParticleName() == "sexaq") {
108  edm::LogVerbatim("CustomPhysics") << " No pmanager implemented for sexaq, only for anti_sexaq";
109  }
110  }
111  }
112  }
113 }
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
static bool s_isgluinoHadron(int pdg)
edm::ParameterSet myConfig
std::string particleDefFilePath
G4ParticleDefinition * GetCloud()
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
void ConstructParticle() override
G4ParticleDefinition * GetSpectator()
std::unique_ptr< CustomParticleFactory > fParticleFactory
static bool s_issbottomHadron(int pdg)
static bool s_isstopHadron(int pdg)
~CustomPhysicsList() override