CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes | Static Private Attributes
CustomPhysicsListSS Class Reference

#include <CustomPhysicsListSS.h>

Inheritance diagram for CustomPhysicsListSS:

Public Member Functions

void ConstructParticle () override
 
void ConstructProcess () override
 
 CustomPhysicsListSS (const std::string &name, const edm::ParameterSet &p, bool useuni=false)
 
 ~CustomPhysicsListSS () override
 

Private Attributes

double dfactor
 
bool fHadronicInteraction
 
std::unique_ptr< CustomParticleFactoryfParticleFactory
 
edm::ParameterSet myConfig
 
std::string particleDefFilePath
 
std::string processDefFilePath
 

Static Private Attributes

static G4ThreadLocal std::unique_ptr< G4Decay > fDecayProcess
 
static G4ThreadLocal std::unique_ptr< G4ProcessHelpermyHelper
 

Detailed Description

Definition at line 12 of file CustomPhysicsListSS.h.

Constructor & Destructor Documentation

CustomPhysicsListSS::CustomPhysicsListSS ( const std::string &  name,
const edm::ParameterSet p,
bool  useuni = false 
)

Definition at line 26 of file CustomPhysicsListSS.cc.

References dfactor, fDecayProcess, fHadronicInteraction, fParticleFactory, edm::ParameterSet::getParameter(), myConfig, myHelper, AlCaHLTBitMon_ParallelJobs::p, and particleDefFilePath.

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 }
std::string particleDefFilePath
T getParameter(std::string const &) const
std::unique_ptr< CustomParticleFactory > fParticleFactory
static G4ThreadLocal std::unique_ptr< G4Decay > fDecayProcess
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
edm::ParameterSet myConfig
CustomPhysicsListSS::~CustomPhysicsListSS ( )
override

Definition at line 50 of file CustomPhysicsListSS.cc.

50  {
51 }

Member Function Documentation

void CustomPhysicsListSS::ConstructParticle ( )
override

Definition at line 53 of file CustomPhysicsListSS.cc.

References fParticleFactory, and particleDefFilePath.

53  {
54  edm::LogInfo("SimG4CoreCustomPhysicsSS")
55  << "===== CustomPhysicsList::ConstructParticle ";
56  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
57 }
std::string particleDefFilePath
std::unique_ptr< CustomParticleFactory > fParticleFactory
void CustomPhysicsListSS::ConstructProcess ( )
override

Definition at line 59 of file CustomPhysicsListSS.cc.

References SimDataFormats::CaloAnalysis::cp, dfactor, fDecayProcess, fHadronicInteraction, fParticleFactory, CustomParticle::GetCloud(), CustomParticle::GetSpectator(), GeV, myConfig, myHelper, and CustomPDGParser::s_isRHadron().

59  {
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 }
const double GeV
Definition: MathUtil.h:16
std::unique_ptr< CustomParticleFactory > fParticleFactory
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)

Member Data Documentation

double CustomPhysicsListSS::dfactor
private

Definition at line 35 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

G4ThreadLocal std::unique_ptr< G4Decay > CustomPhysicsListSS::fDecayProcess
staticprivate

Definition at line 24 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

bool CustomPhysicsListSS::fHadronicInteraction
private

Definition at line 29 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

std::unique_ptr<CustomParticleFactory> CustomPhysicsListSS::fParticleFactory
private

Definition at line 27 of file CustomPhysicsListSS.h.

Referenced by ConstructParticle(), ConstructProcess(), and CustomPhysicsListSS().

edm::ParameterSet CustomPhysicsListSS::myConfig
private

Definition at line 31 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

G4ThreadLocal std::unique_ptr< G4ProcessHelper > CustomPhysicsListSS::myHelper
staticprivate

Definition at line 25 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

std::string CustomPhysicsListSS::particleDefFilePath
private

Definition at line 33 of file CustomPhysicsListSS.h.

Referenced by ConstructParticle(), and CustomPhysicsListSS().

std::string CustomPhysicsListSS::processDefFilePath
private

Definition at line 34 of file CustomPhysicsListSS.h.