Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:30:20

0001 #include <memory>
0002 
0003 #include "SimG4Core/CustomPhysics/interface/CustomPhysicsListSS.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/G4ProcessHelper.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 "G4hMultipleScattering.hh"
0015 #include "G4hIonisation.hh"
0016 #include "G4CoulombScattering.hh"
0017 #include "G4ProcessManager.hh"
0018 
0019 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.h"
0020 #include "SimG4Core/CustomPhysics/interface/CMSDarkPairProductionProcess.h"
0021 
0022 using namespace CLHEP;
0023 
0024 G4ThreadLocal std::unique_ptr<G4ProcessHelper> CustomPhysicsListSS::myHelper;
0025 
0026 CustomPhysicsListSS::CustomPhysicsListSS(const std::string& name, const edm::ParameterSet& p, bool apinew)
0027     : G4VPhysicsConstructor(name) {
0028   myConfig = p;
0029   if (apinew) {
0030     dfactor = p.getParameter<double>("DarkMPFactor");
0031     fHadronicInteraction = p.getParameter<bool>("RhadronPhysics");
0032   } else {
0033     // this is left for backward compatibility
0034     dfactor = p.getParameter<double>("dark_factor");
0035     fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
0036   }
0037   edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
0038   particleDefFilePath = fp.fullPath();
0039   fParticleFactory = std::make_unique<CustomParticleFactory>();
0040   myHelper.reset(nullptr);
0041 
0042   edm::LogVerbatim("SimG4CoreCustomPhysics") << "CustomPhysicsListSS: Path for custom particle definition file: \n"
0043                                              << particleDefFilePath;
0044 }
0045 
0046 CustomPhysicsListSS::~CustomPhysicsListSS() {}
0047 
0048 void CustomPhysicsListSS::ConstructParticle() {
0049   edm::LogVerbatim("SimG4CoreCustomPhysicsSS") << "===== CustomPhysicsList::ConstructParticle ";
0050   fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
0051 }
0052 
0053 void CustomPhysicsListSS::ConstructProcess() {
0054   edm::LogVerbatim("SimG4CoreCustomPhysicsSS") << "CustomPhysicsListSS: adding CustomPhysics processes";
0055 
0056   G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
0057 
0058   for (auto particle : fParticleFactory.get()->getCustomParticles()) {
0059     CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
0060     if (cp) {
0061       G4ProcessManager* pmanager = particle->GetProcessManager();
0062       edm::LogVerbatim("SimG4CoreCustomPhysics")
0063           << "CustomPhysicsListSS: " << particle->GetParticleName() << " PDGcode= " << particle->GetPDGEncoding()
0064           << " Mass= " << particle->GetPDGMass() / GeV << " GeV.";
0065       if (pmanager) {
0066         if (particle->GetPDGCharge() != 0.0) {
0067           ph->RegisterProcess(new G4CoulombScattering, particle);
0068           ph->RegisterProcess(new G4hIonisation, particle);
0069         }
0070         if (cp->GetCloud() && fHadronicInteraction &&
0071             (CustomPDGParser::s_isgluinoHadron(particle->GetPDGEncoding()) ||
0072              (CustomPDGParser::s_isstopHadron(particle->GetPDGEncoding())) ||
0073              (CustomPDGParser::s_issbottomHadron(particle->GetPDGEncoding())))) {
0074           edm::LogVerbatim("SimG4CoreCustomPhysics")
0075               << "CustomPhysicsListSS: " << particle->GetParticleName()
0076               << " CloudMass= " << cp->GetCloud()->GetPDGMass() / GeV
0077               << " GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass() / GeV << " GeV.";
0078 
0079           if (!myHelper.get()) {
0080             myHelper = std::make_unique<G4ProcessHelper>(myConfig, fParticleFactory.get());
0081           }
0082           pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
0083         }
0084         if (particle->GetParticleType() == "darkpho") {
0085           CMSDarkPairProductionProcess* darkGamma = new CMSDarkPairProductionProcess(dfactor);
0086           pmanager->AddDiscreteProcess(darkGamma);
0087         }
0088       }
0089     }
0090   }
0091 }