CMS 3D CMS Logo

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

#include <CustomPhysicsList.h>

Inheritance diagram for CustomPhysicsList:

Public Member Functions

void ConstructParticle () override
 
void ConstructProcess () override
 
 CustomPhysicsList (const std::string &name, const edm::ParameterSet &p)
 
 ~CustomPhysicsList () 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 13 of file CustomPhysicsList.h.

Constructor & Destructor Documentation

CustomPhysicsList::CustomPhysicsList ( const std::string &  name,
const edm::ParameterSet p 
)

Definition at line 25 of file CustomPhysicsList.cc.

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

26  : G4VPhysicsConstructor(name)
27 {
28  myConfig = p;
29  dfactor = p.getParameter<double>("dark_factor");
30  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
31  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
32 
33  particleDefFilePath = fp.fullPath();
34 
36  fDecayProcess.reset(nullptr);
37  myHelper.reset(nullptr);
38 
39  edm::LogInfo("SimG4CoreCustomPhysics")
40  << "CustomPhysicsList: Path for custom particle definition file: \n"
41  <<particleDefFilePath << "\n" << " dark_factor= " << dfactor;
42 }
T getParameter(std::string const &) const
static G4ThreadLocal std::unique_ptr< G4Decay > fDecayProcess
edm::ParameterSet myConfig
std::string particleDefFilePath
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
std::unique_ptr< CustomParticleFactory > fParticleFactory
CustomPhysicsList::~CustomPhysicsList ( )
override

Definition at line 44 of file CustomPhysicsList.cc.

44  {
45 }

Member Function Documentation

void CustomPhysicsList::ConstructParticle ( )
override

Definition at line 47 of file CustomPhysicsList.cc.

References fParticleFactory, and particleDefFilePath.

47  {
48  edm::LogInfo("SimG4CoreCustomPhysics")
49  << "===== CustomPhysicsList::ConstructParticle ";
50  fParticleFactory.get()->loadCustomParticles(particleDefFilePath);
51 }
std::string particleDefFilePath
std::unique_ptr< CustomParticleFactory > fParticleFactory
void CustomPhysicsList::ConstructProcess ( )
override

Definition at line 53 of file CustomPhysicsList.cc.

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

53  {
54 
55  edm::LogInfo("SimG4CoreCustomPhysics")
56  <<"CustomPhysicsList: adding CustomPhysics processes "
57  << "for the list of particles";
58 
59  fDecayProcess.reset(new G4Decay());
60  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
61 
62  for(auto particle : fParticleFactory.get()->GetCustomParticles()) {
63  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
64  if(cp) {
65  G4ProcessManager* pmanager = particle->GetProcessManager();
66  edm::LogInfo("SimG4CoreCustomPhysics")
67  <<"CustomPhysicsList: " << particle->GetParticleName()
68  <<" PDGcode= " << particle->GetPDGEncoding()
69  << " Mass= " << particle->GetPDGMass()/GeV <<" GeV.";
70  if(pmanager) {
71  if(particle->GetPDGCharge() != 0.0) {
72  ph->RegisterProcess(new G4hMultipleScattering, particle);
73  ph->RegisterProcess(new G4hIonisation, particle);
74  }
75  if(fDecayProcess.get()->IsApplicable(*particle)) {
76  ph->RegisterProcess(fDecayProcess.get(), particle);
77  }
78  if(cp->GetCloud() && fHadronicInteraction &&
79  CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
80  edm::LogInfo("SimG4CoreCustomPhysics")
81  <<"CustomPhysicsList: " << particle->GetParticleName()
82  <<" CloudMass= " <<cp->GetCloud()->GetPDGMass()/GeV
83  <<" GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass()/GeV<<" GeV.";
84 
85  if(!myHelper.get()) { myHelper.reset(new G4ProcessHelper(myConfig, fParticleFactory.get())); }
86  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper.get()));
87  }
88  if(particle->GetParticleType()=="darkpho"){
90  pmanager->AddDiscreteProcess(darkGamma);
91  }
92  }
93  }
94  }
95 }
const double GeV
Definition: MathUtil.h:16
static G4ThreadLocal std::unique_ptr< G4Decay > fDecayProcess
edm::ParameterSet myConfig
G4ParticleDefinition * GetCloud()
static G4ThreadLocal std::unique_ptr< G4ProcessHelper > myHelper
G4ParticleDefinition * GetSpectator()
static bool s_isRHadron(int pdg)
std::unique_ptr< CustomParticleFactory > fParticleFactory

Member Data Documentation

double CustomPhysicsList::dfactor
private

Definition at line 35 of file CustomPhysicsList.h.

Referenced by ConstructProcess(), and CustomPhysicsList().

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

Definition at line 24 of file CustomPhysicsList.h.

Referenced by ConstructProcess(), and CustomPhysicsList().

bool CustomPhysicsList::fHadronicInteraction
private

Definition at line 29 of file CustomPhysicsList.h.

Referenced by ConstructProcess(), and CustomPhysicsList().

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

Definition at line 27 of file CustomPhysicsList.h.

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

edm::ParameterSet CustomPhysicsList::myConfig
private

Definition at line 31 of file CustomPhysicsList.h.

Referenced by ConstructProcess(), and CustomPhysicsList().

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

Definition at line 25 of file CustomPhysicsList.h.

Referenced by ConstructProcess(), and CustomPhysicsList().

std::string CustomPhysicsList::particleDefFilePath
private

Definition at line 33 of file CustomPhysicsList.h.

Referenced by ConstructParticle(), and CustomPhysicsList().

std::string CustomPhysicsList::processDefFilePath
private

Definition at line 34 of file CustomPhysicsList.h.