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