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