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"
5 
9 
10 #include "G4Decay.hh"
11 #include "G4hMultipleScattering.hh"
12 #include "G4hIonisation.hh"
13 #include "G4ProcessManager.hh"
14 
15 #include "G4LeptonConstructor.hh"
16 #include "G4MesonConstructor.hh"
17 #include "G4BaryonConstructor.hh"
18 #include "G4ShortLivedConstructor.hh"
19 #include "G4IonConstructor.hh"
20 
21 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
22 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
23 
24 using namespace CLHEP;
25 
26 
28  : G4VPhysicsConstructor(name)
29 {
30  myConfig = p;
31  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
33  edm::LogInfo("CustomPhysics")<<"Path for custom particle definition file: "
35  myHelper = 0;
36 }
37 
39  delete myHelper;
40 }
41 
44 }
45 
48 }
49 
51 
52  edm::LogInfo("CustomPhysics") << " CustomPhysicsList: adding CustomPhysics processes "
53  << "for the list of particles: \n";
54  aParticleIterator->reset();
55 
56  while((*aParticleIterator)()) {
57  G4ParticleDefinition* particle = aParticleIterator->value();
59  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
60  edm::LogInfo("CustomPhysics") << particle->GetParticleName()
61  <<" PDGcode= "<<particle->GetPDGEncoding()
62  << " Mass= "
63  <<particle->GetPDGMass()/GeV <<" GeV.";
64  if(cp->GetCloud()!=0) {
65  edm::LogInfo("CustomPhysics") << particle->GetParticleName()
66  <<" CloudMass= "
67  <<cp->GetCloud()->GetPDGMass()/GeV
68  <<" GeV; SpectatorMass= "
69  << cp->GetSpectator()->GetPDGMass()/GeV
70  <<" GeV.";
71  }
72  G4ProcessManager* pmanager = particle->GetProcessManager();
73  if(pmanager) {
74  if(particle->GetPDGCharge() != 0.0) {
75  pmanager->AddProcess(new G4hMultipleScattering,-1, 1, 1);
76  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
77  }
78  if(cp != 0) {
79  if(particle->GetParticleType()=="rhadron" ||
80  particle->GetParticleType()=="mesonino" ||
81  particle->GetParticleType() == "sbaryon"){
82  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
83  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper));
84  }
85  }
86  }
87  else LogDebug("CustomPhysics") << " No pmanager";
88  }
89  }
90 }
91 
92 
93 void CustomPhysicsList::setupRHadronPhycis(G4ParticleDefinition* particle)
94 {
95  // LogDebug("CustomPhysics")<<"Configuring rHadron: "
96  // <<cp->
97 
98  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
99  if(cp->GetCloud()!=0) {
100  LogDebug("CustomPhysics")
101  <<"Cloud mass is "
102  <<cp->GetCloud()->GetPDGMass()/GeV
103  <<" GeV. Spectator mass is "
104  <<static_cast<CustomParticle*>(particle)->GetSpectator()->GetPDGMass()/GeV
105  <<" GeV.";
106  }
107  G4ProcessManager* pmanager = particle->GetProcessManager();
108  if(pmanager){
109  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
110  if(particle->GetPDGCharge()/eplus != 0){
111  pmanager->AddProcess(new G4hMultipleScattering,-1, 1, 1);
112  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
113  }
114  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper)); //GHEISHA
115  }
116  else LogDebug("CustomPhysics") << " No pmanager";
117 }
118 
119 void CustomPhysicsList::setupSUSYPhycis(G4ParticleDefinition* particle)
120 {
121  G4ProcessManager* pmanager = particle->GetProcessManager();
122  if(pmanager){
123  if(particle->GetPDGCharge()/eplus != 0){
124  pmanager->AddProcess(new G4hMultipleScattering,-1, 1, 1);
125  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
126  }
127  pmanager->AddProcess(new G4Decay, 1, -1, 3);
128  }
129  else LogDebug("CustomPhysics") << " No pmanager";
130 }
#define LogDebug(id)
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
edm::ParameterSet myConfig
std::string particleDefFilePath
G4ParticleDefinition * GetCloud()
void setupRHadronPhycis(G4ParticleDefinition *particle)
static void loadCustomParticles(const std::string &filePath)
virtual ~CustomPhysicsList()
void setupSUSYPhycis(G4ParticleDefinition *particle)
G4ParticleDefinition * GetSpectator()
G4ProcessHelper * myHelper
CustomPhysicsList(std::string name, const edm::ParameterSet &p)
static bool isCustomParticle(G4ParticleDefinition *particle)
std::string fullPath() const
Definition: FileInPath.cc:165
virtual void ConstructParticle()
virtual void ConstructProcess()