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"
6 
10 
11 #include "G4Decay.hh"
12 #include "G4hMultipleScattering.hh"
13 #include "G4hIonisation.hh"
14 #include "G4ProcessManager.hh"
15 
16 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
17 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
18 #include "SimG4Core/CustomPhysics/interface/CMSDarkPairProductionProcess.hh"
19 
20 using namespace CLHEP;
21 
22 G4ThreadLocal G4Decay* CustomPhysicsList::fDecayProcess = nullptr;
23 G4ThreadLocal G4ProcessHelper* CustomPhysicsList::myHelper = nullptr;
24 
26  : G4VPhysicsConstructor(name)
27 {
28  myConfig = p;
29  dfactor = p.getParameter<double>("dark_factor");
30  edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
31  fHadronicInteraction = p.getParameter<bool>("rhadronPhysics");
32 
33  particleDefFilePath = fp.fullPath();
34  edm::LogInfo("SimG4CoreCustomPhysics")
35  << "CustomPhysicsList: Path for custom particle definition file: \n"
36  <<particleDefFilePath << "\n" << " dark_factor= " << dfactor;
37 }
38 
40 }
41 
43  G4cout << "===== CustomPhysicsList::ConstructParticle " << this << G4endl;
45 }
46 
48 
49  edm::LogInfo("SimG4CoreCustomPhysics")
50  <<"CustomPhysicsList: adding CustomPhysics processes "
51  << "for the list of particles";
52 
53  fDecayProcess = new G4Decay();
54 
55  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
56 
57  aParticleIterator->reset();
58  G4ParticleDefinition* particle;
59 
60  while((*aParticleIterator)()) {
61  particle = aParticleIterator->value();
63  CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
64  G4ProcessManager* pmanager = particle->GetProcessManager();
65  edm::LogInfo("SimG4CoreCustomPhysics")
66  <<"CustomPhysicsList: " << particle->GetParticleName()
67  <<" PDGcode= " << particle->GetPDGEncoding()
68  << " Mass= " << particle->GetPDGMass()/GeV <<" GeV.";
69  if(cp && pmanager) {
70  if(particle->GetPDGCharge() != 0.0) {
71  ph->RegisterProcess(new G4hMultipleScattering, particle);
72  ph->RegisterProcess(new G4hIonisation, particle);
73  }
74  if(fDecayProcess->IsApplicable(*particle)) {
75  ph->RegisterProcess(fDecayProcess, particle);
76  }
77  if(cp->GetCloud() && fHadronicInteraction &&
78  CustomPDGParser::s_isRHadron(particle->GetPDGEncoding())) {
79  edm::LogInfo("SimG4CoreCustomPhysics")
80  <<"CustomPhysicsList: " << particle->GetParticleName()
81  <<" CloudMass= " <<cp->GetCloud()->GetPDGMass()/GeV
82  <<" GeV; SpectatorMass= " << cp->GetSpectator()->GetPDGMass()/GeV<<" GeV.";
83 
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  }
94 }
T getParameter(std::string const &) const
const double GeV
Definition: MathUtil.h:16
static G4ThreadLocal G4Decay * fDecayProcess
edm::ParameterSet myConfig
std::string particleDefFilePath
G4ParticleDefinition * GetCloud()
static void loadCustomParticles(const std::string &filePath)
virtual ~CustomPhysicsList()
static G4ThreadLocal G4ProcessHelper * myHelper
G4ParticleDefinition * GetSpectator()
static bool s_isRHadron(int pdg)
CustomPhysicsList(std::string name, const edm::ParameterSet &p)
static bool isCustomParticle(G4ParticleDefinition *particle)
virtual void ConstructParticle()
virtual void ConstructProcess()