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 #include "SimG4Core/CustomPhysics/interface/CMSDarkPairProductionProcess.hh"
24 
25 using namespace CLHEP;
26 
27 
29  : G4VPhysicsConstructor(name)
30 {
31  myConfig = p;
32  dfactor = p.getParameter<double>("dark_factor");
33  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
34  particleDefFilePath = fp.fullPath();
35  edm::LogInfo("CustomPhysics")<<"Path for custom particle definition file: "
36  <<particleDefFilePath<< "\n" << " dark_factor= " << dfactor;
37  myHelper = 0;
38 }
39 
41  delete myHelper;
42 }
43 
46 }
47 
50 }
51 
53 
54  edm::LogInfo("CustomPhysics") << " CustomPhysicsList: adding CustomPhysics processes "
55  << "for the list of particles: \n";
56  aParticleIterator->reset();
57 
58  while((*aParticleIterator)()) {
59  G4ParticleDefinition* particle = aParticleIterator->value();
61  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
62  edm::LogInfo("CustomPhysics") << particle->GetParticleName()
63  <<" PDGcode= "<<particle->GetPDGEncoding()
64  << " Mass= "
65  <<particle->GetPDGMass()/GeV <<" GeV.";
66  if(cp->GetCloud()!=0) {
67  edm::LogInfo("CustomPhysics") << particle->GetParticleName()
68  <<" CloudMass= "
69  <<cp->GetCloud()->GetPDGMass()/GeV
70  <<" GeV; SpectatorMass= "
71  << cp->GetSpectator()->GetPDGMass()/GeV
72  <<" GeV.";
73  }
74  G4ProcessManager* pmanager = particle->GetProcessManager();
75  if(pmanager) {
76  if(particle->GetPDGCharge() != 0.0) {
77  pmanager->AddProcess(new G4hMultipleScattering,-1, 1, 1);
78  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
79  }
80  if(cp != 0) {
81  if(particle->GetParticleType()=="rhadron" ||
82  particle->GetParticleType()=="mesonino" ||
83  particle->GetParticleType() == "sbaryon"){
84  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
85  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper));
86  }
87  if(particle->GetParticleType()=="darkpho"){
88  CMSDarkPairProductionProcess * darkGamma = new CMSDarkPairProductionProcess(dfactor);
89  pmanager->AddDiscreteProcess(darkGamma);
90  }
91  }
92  }
93  else LogDebug("CustomPhysics") << " No pmanager";
94  }
95  }
96 }
97 
98 
99 void CustomPhysicsList::setupRHadronPhycis(G4ParticleDefinition* particle)
100 {
101  // LogDebug("CustomPhysics")<<"Configuring rHadron: "
102  // <<cp->
103 
104  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
105  if(cp->GetCloud()!=0) {
106  LogDebug("CustomPhysics")
107  <<"Cloud mass is "
108  <<cp->GetCloud()->GetPDGMass()/GeV
109  <<" GeV. Spectator mass is "
110  <<static_cast<CustomParticle*>(particle)->GetSpectator()->GetPDGMass()/GeV
111  <<" GeV.";
112  }
113  G4ProcessManager* pmanager = particle->GetProcessManager();
114  if(pmanager){
115  if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
116  if(particle->GetPDGCharge()/eplus != 0){
117  pmanager->AddProcess(new G4hMultipleScattering,-1, 1, 1);
118  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
119  }
120  pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper)); //GHEISHA
121  }
122  else LogDebug("CustomPhysics") << " No pmanager";
123 }
124 
125 void CustomPhysicsList::setupSUSYPhycis(G4ParticleDefinition* particle)
126 {
127  G4ProcessManager* pmanager = particle->GetProcessManager();
128  if(pmanager){
129  if(particle->GetPDGCharge()/eplus != 0){
130  pmanager->AddProcess(new G4hMultipleScattering,-1, 1, 1);
131  pmanager->AddProcess(new G4hIonisation, -1, 2, 2);
132  }
133  pmanager->AddProcess(new G4Decay, 1, -1, 3);
134  }
135  else LogDebug("CustomPhysics") << " No pmanager";
136 }
#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)
virtual void ConstructParticle()
virtual void ConstructProcess()