CMS 3D CMS Logo

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