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 (std::string name, const edm::ParameterSet &p)
 
 ~CustomPhysicsListSS () override
 

Private Attributes

bool fHadronicInteraction
 
edm::ParameterSet myConfig
 
std::string particleDefFilePath
 
std::string processDefFilePath
 

Static Private Attributes

static G4ThreadLocal G4Decay * fDecayProcess = 0
 
static G4ThreadLocal G4ProcessHelper * myHelper = 0
 

Detailed Description

Definition at line 15 of file CustomPhysicsListSS.h.

Constructor & Destructor Documentation

CustomPhysicsListSS::CustomPhysicsListSS ( std::string  name,
const edm::ParameterSet p 
)

Definition at line 25 of file CustomPhysicsListSS.cc.

References fHadronicInteraction, edm::FileInPath::fullPath(), edm::ParameterSet::getParameter(), myConfig, myHelper, AlCaHLTBitMon_ParallelJobs::p, and particleDefFilePath.

26  : G4VPhysicsConstructor(name)
27 {
28  myConfig = p;
29  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
30  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
32  edm::LogInfo("SimG4CoreCustomPhysics")
33  <<"CustomPhysicsListSS: Path for custom particle definition file: \n"
35  myHelper = nullptr;
36 }
std::string particleDefFilePath
T getParameter(std::string const &) const
static G4ThreadLocal G4ProcessHelper * myHelper
edm::ParameterSet myConfig
std::string fullPath() const
Definition: FileInPath.cc:184
CustomPhysicsListSS::~CustomPhysicsListSS ( )
override

Definition at line 38 of file CustomPhysicsListSS.cc.

38  {
39 }

Member Function Documentation

void CustomPhysicsListSS::ConstructParticle ( )
override

Definition at line 41 of file CustomPhysicsListSS.cc.

References CustomParticleFactory::loadCustomParticles(), and particleDefFilePath.

41  {
43 }
std::string particleDefFilePath
static void loadCustomParticles(const std::string &filePath)
void CustomPhysicsListSS::ConstructProcess ( )
override

Definition at line 45 of file CustomPhysicsListSS.cc.

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

45  {
46 
47  edm::LogInfo("SimG4CoreCustomPhysics")
48  <<"CustomPhysicsListSS: adding CustomPhysics processes";
49 
50  fDecayProcess = new G4Decay();
51  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
52 
53  for(auto particle : CustomParticleFactory::GetCustomParticles()) {
54 
55  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
56  if(cp) {
57  G4ProcessManager* pmanager = particle->GetProcessManager();
58  edm::LogInfo("SimG4CoreCustomPhysics")
59  <<"CustomPhysicsListSS: " << particle->GetParticleName()
60  <<" PDGcode= " << particle->GetPDGEncoding()
61  << " Mass= " << particle->GetPDGMass()/GeV <<" GeV.";
62  if(pmanager) {
63  if(particle->GetPDGCharge() != 0.0) {
64  ph->RegisterProcess(new G4hMultipleScattering, particle);
65  ph->RegisterProcess(new G4hIonisation, particle);
66  }
67  if(fDecayProcess->IsApplicable(*particle)) {
68  ph->RegisterProcess(fDecayProcess, particle);
69  }
70  if(cp->GetCloud() && fHadronicInteraction &&
71  CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
72  edm::LogInfo("SimG4CoreCustomPhysics")
73  <<"CustomPhysicsListSS: " << particle->GetParticleName()
74  <<" CloudMass= " <<cp->GetCloud()->GetPDGMass()/GeV
75  <<" GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass()/GeV <<" GeV.";
76 
77  if(!myHelper) { myHelper = new G4ProcessHelper(myConfig); }
78  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper));
79  }
80  }
81  }
82  }
83 }
const double GeV
Definition: MathUtil.h:16
G4ParticleDefinition * GetCloud()
static G4ThreadLocal G4ProcessHelper * myHelper
static G4ThreadLocal G4Decay * fDecayProcess
G4ParticleDefinition * GetSpectator()
edm::ParameterSet myConfig
static bool s_isRHadron(int pdg)
static const std::vector< G4ParticleDefinition * > & GetCustomParticles()

Member Data Documentation

G4ThreadLocal G4Decay * CustomPhysicsListSS::fDecayProcess = 0
staticprivate

Definition at line 26 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess().

bool CustomPhysicsListSS::fHadronicInteraction
private

Definition at line 29 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

edm::ParameterSet CustomPhysicsListSS::myConfig
private

Definition at line 31 of file CustomPhysicsListSS.h.

Referenced by ConstructProcess(), and CustomPhysicsListSS().

G4ThreadLocal G4ProcessHelper * CustomPhysicsListSS::myHelper = 0
staticprivate

Definition at line 27 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.