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 using namespace CLHEP;
24 
25 
26 CustomPhysicsList::CustomPhysicsList(std::string name, const edm::ParameterSet & p) : G4VPhysicsConstructor(name) {
27 
28  myConfig = p;
29  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
31  edm::LogInfo("CustomPhysics")<<"Path for custom particle definition file: "<<particleDefFilePath;
32  myHelper = 0;
33 
34  }
35 
37  delete myHelper;
38 }
39 
42 }
43 
46 }
47 
49  LogDebug("CustomPhysics") << " CustomPhysics: adding CustomPhysics processes";
50  theParticleIterator->reset();
51 
52  while((*theParticleIterator)()) {
53  int i = 0;
54  G4ParticleDefinition* particle = theParticleIterator->value();
55  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
57  LogDebug("CustomPhysics") << particle->GetParticleName()<<", "<<particle->GetPDGEncoding()
58  << " is Custom. Mass is "
59  <<particle->GetPDGMass()/GeV <<" GeV.";
60  if(cp->GetCloud()!=0) {
61  LogDebug("CustomPhysics")<<"Cloud mass is "
62  <<cp->GetCloud()->GetPDGMass()/GeV
63  <<" GeV. Spectator mass is "
64  <<static_cast<CustomParticle*>(particle)->GetSpectator()->GetPDGMass()/GeV
65  <<" GeV.";
66  }
67  G4ProcessManager* pmanager = particle->GetProcessManager();
68  if(pmanager) {
69  if(cp!=0) {
70  if(particle->GetParticleType()=="rhadron" || particle->GetParticleType()=="mesonino" || particle->GetParticleType() == "sbaryon"){
71  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
72  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper)); //GHEISHA
73  }
74  }
75  if(particle->GetPDGCharge()/eplus != 0) {
76  pmanager->AddProcess(new G4hMultipleScattering,-1, 1,i+1);
77  pmanager->AddProcess(new G4hIonisation, -1, 2,i+2);
78  }
79  }
80  else LogDebug("CustomPhysics") << " No pmanager";
81  }
82  }
83 }
84 
85 
86 void CustomPhysicsList::setupRHadronPhycis(G4ParticleDefinition* particle){
87 
88  // LogDebug("CustomPhysics")<<"Configuring rHadron: "
89  // <<cp->
90 
91  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
92  if(cp->GetCloud()!=0)
93  LogDebug("CustomPhysics")<<"Cloud mass is "
94  <<cp->GetCloud()->GetPDGMass()/GeV
95  <<" GeV. Spectator mass is "
96  <<static_cast<CustomParticle*>(particle)->GetSpectator()->GetPDGMass()/GeV
97  <<" GeV.";
98 
99  G4ProcessManager* pmanager = particle->GetProcessManager();
100  if(pmanager){
101  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
102  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper)); //GHEISHA
103  if(particle->GetPDGCharge()/eplus != 0){
104  pmanager->AddProcess(new G4hMultipleScattering,-1, 1,1);
105  pmanager->AddProcess(new G4hIonisation, -1, 2,2);
106  }
107  }
108  else LogDebug("CustomPhysics") << " No pmanager";
109 }
110 
111 
112 void CustomPhysicsList::setupSUSYPhycis(G4ParticleDefinition* particle){
113 
114 // CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
115  G4ProcessManager* pmanager = particle->GetProcessManager();
116  if(pmanager){
117  pmanager->AddProcess(new G4Decay,1, 1,1);
118  if(particle->GetPDGCharge()/eplus != 0){
119  pmanager->AddProcess(new G4hMultipleScattering,-1, 2,2);
120  pmanager->AddProcess(new G4hIonisation, -1, 3,3);
121  }
122  }
123  else LogDebug("CustomPhysics") << " No pmanager";
124 }
#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()