Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-07-09 05:00:33

0001 #include <memory>
0002 
0003 #include "SimG4Core/CustomPhysics/interface/CustomPhysicsList.h"
0004 #include "SimG4Core/CustomPhysics/interface/CustomParticleFactory.h"
0005 #include "SimG4Core/CustomPhysics/interface/CustomParticle.h"
0006 #include "SimG4Core/CustomPhysics/interface/DummyChargeFlipProcess.h"
0007 #include "SimG4Core/CustomPhysics/interface/CustomProcessHelper.h"
0008 #include "SimG4Core/CustomPhysics/interface/CustomPDGParser.h"
0009 
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/ParameterSet/interface/FileInPath.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 
0014 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.h"
0015 #include "SimG4Core/CustomPhysics/interface/CMSDarkPairProductionProcess.h"
0016 #include "SimG4Core/CustomPhysics/interface/CMSQGSPSIMPBuilder.h"
0017 #include "SimG4Core/CustomPhysics/interface/CMSSIMPInelasticProcess.h"
0018 
0019 #include "SimG4Core/CustomPhysics/interface/CMSSQLoopProcess.h"
0020 #include "SimG4Core/CustomPhysics/interface/CMSSQLoopProcessDiscr.h"
0021 #include "SimG4Core/CustomPhysics/interface/CMSSQNeutronAnnih.h"
0022 #include "SimG4Core/CustomPhysics/interface/CMSSQInelasticCrossSection.h"
0023 
0024 #include "G4hMultipleScattering.hh"
0025 #include "G4hIonisation.hh"
0026 #include "G4ProcessManager.hh"
0027 #include "G4HadronicProcess.hh"
0028 #include "G4AutoDelete.hh"
0029 
0030 using namespace CLHEP;
0031 
0032 G4ThreadLocal CustomProcessHelper* CustomPhysicsList::myHelper = nullptr;
0033 
0034 CustomPhysicsList::CustomPhysicsList(const std::string& name, const edm::ParameterSet& p, bool apinew)
0035     : G4VPhysicsConstructor(name) {
0036   myConfig = p;
0037   if (apinew) {
0038     dfactor = p.getParameter<double>("DarkMPFactor");
0039     fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
0040   } else {
0041     // this is left for backward compatibility
0042     dfactor = p.getParameter<double>("dark_factor");
0043     fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
0044   }
0045   edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
0046   particleDefFilePath = fp.fullPath();
0047   fParticleFactory = std::make_unique<CustomParticleFactory>();
0048 
0049   edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsList: Path for custom particle definition file: \n"
0050                                              << particleDefFilePath << "\n"
0051                                              << "      dark_factor= " << dfactor;
0052 }
0053 
0054 void CustomPhysicsList::ConstructParticle() {
0055   edm::LogVerbatim("SimG4CoreCustomPhysics") << "===== CustomPhysicsList::ConstructParticle ";
0056   fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
0057 }
0058 
0059 void CustomPhysicsList::ConstructProcess() {
0060   edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsList: adding CustomPhysics processes "
0061                                              << "for the list of particles";
0062 
0063   G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
0064 
0065   for (auto particle : fParticleFactory.get()->getCustomParticles()) {
0066     if (particle->GetParticleType() == "simp") {
0067       G4ProcessManager* pmanager = particle->GetProcessManager();
0068       if (pmanager) {
0069         CMSSIMPInelasticProcess* simpInelPr = new CMSSIMPInelasticProcess();
0070         CMSQGSPSIMPBuilder* theQGSPSIMPB = new CMSQGSPSIMPBuilder();
0071         theQGSPSIMPB->Build(simpInelPr);
0072         pmanager->AddDiscreteProcess(simpInelPr);
0073       } else
0074         edm::LogVerbatim("CustomPhysics") << "   No pmanager";
0075     }
0076 
0077     CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
0078     if (cp) {
0079       G4ProcessManager* pmanager = particle->GetProcessManager();
0080       edm::LogVerbatim("SimG4CoreCustomPhysics")
0081           << "CustomPhysicsList: " << particle->GetParticleName() << "  PDGcode= " << particle->GetPDGEncoding()
0082           << "  Mass= " << particle->GetPDGMass() / GeV << " GeV.";
0083       if (pmanager) {
0084         if (particle->GetPDGCharge() != 0.0) {
0085           ph->RegisterProcess(new G4hMultipleScattering, particle);
0086           ph->RegisterProcess(new G4hIonisation, particle);
0087         }
0088         if (cp->GetCloud() && fHadronicInteraction &&
0089             (CustomPDGParser::s_isgluinoHadron(particle->GetPDGEncoding()) ||
0090              (CustomPDGParser::s_isstopHadron(particle->GetPDGEncoding())) ||
0091              (CustomPDGParser::s_issbottomHadron(particle->GetPDGEncoding())))) {
0092           edm::LogVerbatim("SimG4CoreCustomPhysics")
0093               << "CustomPhysicsList: " << particle->GetParticleName()
0094               << " CloudMass= " << cp->GetCloud()->GetPDGMass() / GeV
0095               << " GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass() / GeV << " GeV.";
0096 
0097           if (nullptr == myHelper) {
0098             myHelper = new CustomProcessHelper(myConfig, fParticleFactory.get());
0099             G4AutoDelete::Register(myHelper);
0100           }
0101           pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper));
0102         }
0103         if (particle->GetParticleType() == "darkpho") {
0104           CMSDarkPairProductionProcess* darkGamma = new CMSDarkPairProductionProcess(dfactor);
0105           pmanager->AddDiscreteProcess(darkGamma);
0106         }
0107         if (particle->GetParticleName() == "anti_sexaq") {
0108           // here the different sexaquark interactions get defined
0109           G4HadronicProcess* sqInelPr = new G4HadronicProcess();
0110           CMSSQNeutronAnnih* sqModel = new CMSSQNeutronAnnih(particle->GetPDGMass() / GeV);
0111           sqInelPr->RegisterMe(sqModel);
0112           CMSSQInelasticCrossSection* sqInelXS = new CMSSQInelasticCrossSection(particle->GetPDGMass() / GeV);
0113           sqInelPr->AddDataSet(sqInelXS);
0114           pmanager->AddDiscreteProcess(sqInelPr);
0115           // add also the looping needed to simulate flat interaction probability
0116           CMSSQLoopProcess* sqLoopPr = new CMSSQLoopProcess();
0117           pmanager->AddContinuousProcess(sqLoopPr);
0118           CMSSQLoopProcessDiscr* sqLoopPrDiscr = new CMSSQLoopProcessDiscr(particle->GetPDGMass() / GeV);
0119           pmanager->AddDiscreteProcess(sqLoopPrDiscr);
0120         } else if (particle->GetParticleName() == "sexaq") {
0121           edm::LogVerbatim("CustomPhysics") << "   No pmanager implemented for sexaq, only for anti_sexaq";
0122         }
0123       }
0124     }
0125   }
0126 }