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