CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CustomPhysicsListSS.cc
Go to the documentation of this file.
4 #include "SimG4Core/CustomPhysics/interface/G4ProcessHelper.hh"
6 
10 
11 #include "G4Decay.hh"
12 #include "G4hMultipleScattering.hh"
13 #include "G4hIonisation.hh"
14 #include "G4CoulombScattering.hh"
15 #include "G4ProcessManager.hh"
16 
17 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
18 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
19 
20 using namespace CLHEP;
21 
22 G4ThreadLocal G4Decay* CustomPhysicsListSS::fDecayProcess = 0;
23 G4ThreadLocal G4ProcessHelper* CustomPhysicsListSS::myHelper = 0;
24 G4ThreadLocal bool CustomPhysicsListSS::fInitialized = false;
25 
27  : G4VPhysicsConstructor(name)
28 {
29  myConfig = p;
30  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
31  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
33  edm::LogInfo("SimG4CoreCustomPhysics")
34  <<"CustomPhysicsListSS: Path for custom particle definition file: \n"
36  myHelper = 0;
37 }
38 
40 }
41 
44 }
45 
47 
48  if(fInitialized) { return; }
49  fInitialized = true;
50 
51  edm::LogInfo("SimG4CoreCustomPhysics")
52  <<"CustomPhysicsListSS: adding CustomPhysics processes";
53 
54  fDecayProcess = new G4Decay();
55 
56  aParticleIterator->reset();
57 
58  while((*aParticleIterator)()) {
59  G4ParticleDefinition* particle = aParticleIterator->value();
61  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
62  G4ProcessManager* pmanager = particle->GetProcessManager();
63  edm::LogInfo("SimG4CoreCustomPhysics")
64  <<"CustomPhysicsListSS: " << particle->GetParticleName()
65  <<" PDGcode= " << particle->GetPDGEncoding()
66  << " Mass= " << particle->GetPDGMass()/GeV <<" GeV.";
67  if(cp && pmanager) {
68  if(particle->GetPDGCharge()/eplus != 0) {
69  pmanager->AddProcess(new G4CoulombScattering, -1,-1, 1);
70  pmanager->AddProcess(new G4hIonisation, -1, 1, 2);
71  }
72  if(fDecayProcess->IsApplicable(*particle)) {
73  pmanager->AddProcess(new G4Decay, 0, -1, 3);
74  }
75  if(cp->GetCloud() && fHadronicInteraction &&
76  CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
77  edm::LogInfo("SimG4CoreCustomPhysics")
78  <<"CustomPhysicsListSS: " << particle->GetParticleName()
79  <<" CloudMass= " <<cp->GetCloud()->GetPDGMass()/GeV
80  <<" GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass()/GeV <<" GeV.";
81 
82  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
83  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper));
84  }
85  }
86  }
87  }
88 }
std::string particleDefFilePath
static G4ThreadLocal bool fInitialized
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
G4ParticleDefinition * GetCloud()
virtual void ConstructParticle()
CustomPhysicsListSS(std::string name, const edm::ParameterSet &p)
static G4ThreadLocal G4ProcessHelper * myHelper
static void loadCustomParticles(const std::string &filePath)
static G4ThreadLocal G4Decay * fDecayProcess
G4ParticleDefinition * GetSpectator()
edm::ParameterSet myConfig
static bool s_isRHadron(int pdg)
static bool isCustomParticle(G4ParticleDefinition *particle)
std::string fullPath() const
Definition: FileInPath.cc:165
virtual void ConstructProcess()