CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
CMSMonopolePhysics.cc
Go to the documentation of this file.
4 
5 #include "G4ParticleDefinition.hh"
6 #include "G4ProcessManager.hh"
7 
8 #include "G4StepLimiter.hh"
9 #include "G4hMultipleScattering.hh"
10 #include "G4hIonisation.hh"
11 #include "G4mplIonisation.hh"
12 #include "G4mplIonisationWithDeltaModel.hh"
13 
14 #include "CLHEP/Units/GlobalSystemOfUnits.h"
15 
17  : G4VPhysicsConstructor("Monopole Physics") {
18  verbose = p.getUntrackedParameter<int>("Verbosity", 0);
19  magCharge = p.getUntrackedParameter<int>("MonopoleCharge", 1);
20  deltaRay = p.getUntrackedParameter<bool>("MonopoleDeltaRay", true);
21  multiSc = p.getUntrackedParameter<bool>("MonopoleMultiScatter", false);
22  transport = p.getUntrackedParameter<bool>("MonopoleTransport", true);
23  double mass = p.getUntrackedParameter<double>("MonopoleMass", 200);
24  if (pdt && mass > 0.0) {
25  int ii = 0;
26  for (HepPDT::ParticleDataTable::const_iterator p = pdt->begin(); p != pdt->end(); ++p, ++ii) {
27  HepPDT::ParticleData particle = (p->second);
28  std::string particleName = (particle.name()).substr(0, 8);
29  if (strcmp(particleName.c_str(), "Monopole") == 0) {
30  names.push_back(particle.name());
31  masses.push_back(mass * CLHEP::GeV);
32  elCharges.push_back((int)(particle.charge()));
33  pdgEncodings.push_back(particle.pid());
34  monopoles.push_back(nullptr);
35  if (verbose > 0)
36  G4cout << "CMSMonopolePhysics: Monopole[" << ii << "] " << particleName << " Mass " << particle.mass()
37  << " GeV, Magnetic Charge " << magCharge << ", Electric Charge " << particle.charge() << G4endl;
38  } else if (strcmp(particleName.c_str(), "AntiMono") == 0) {
39  names.push_back(particle.name());
40  masses.push_back(mass * CLHEP::GeV);
41  elCharges.push_back((int)(particle.charge()));
42  pdgEncodings.push_back(particle.pid());
43  monopoles.push_back(nullptr);
44  if (verbose > 0)
45  G4cout << "CMSMonopolePhysics: Monopole[" << ii << "] " << particleName << " Mass " << particle.mass()
46  << " GeV, Magnetic Charge " << magCharge << ", Electric Charge " << particle.charge() << G4endl;
47  }
48  }
49  }
50  if (verbose > 0)
51  G4cout << "CMSMonopolePhysics has " << names.size() << " monopole candidates and delta Ray option " << deltaRay
52  << G4endl;
53 }
54 
56 
58  for (unsigned int ii = 0; ii < names.size(); ++ii) {
59  // monopoles are created once in the master thread
60  if (!monopoles[ii]) {
61  G4int mc = (pdgEncodings[ii] >= 0) ? magCharge : -magCharge;
62  Monopole* mpl = new Monopole(names[ii], pdgEncodings[ii], masses[ii], mc, elCharges[ii]);
63  monopoles[ii] = mpl;
64  if (verbose > 0)
65  G4cout << "Create Monopole " << names[ii] << " of mass " << masses[ii] / CLHEP::GeV << " GeV, magnetic charge "
66  << mc << ", electric charge " << elCharges[ii] << " and PDG encoding " << pdgEncodings[ii] << " at "
67  << monopoles[ii] << G4endl;
68  }
69  }
70 }
71 
73  // Add standard EM Processes
74  if (verbose > 0) {
75  G4cout << "### CMSMonopolePhysics ConstructProcess()" << G4endl;
76  }
77  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
78 
79  for (unsigned int ii = 0; ii < monopoles.size(); ++ii) {
80  if (monopoles[ii]) {
81  Monopole* mpl = monopoles[ii];
82  G4ProcessManager* pmanager = mpl->GetProcessManager();
83  if (!pmanager) {
84  std::ostringstream o;
85  o << "Monopole without a Process Manager";
86  throw edm::Exception(edm::errors::Configuration, o.str().c_str());
87  return;
88  }
89 
90  G4double magn = mpl->MagneticCharge();
91  G4double mass = mpl->GetPDGMass();
92  if (verbose > 1) {
93  G4cout << "### CMSMonopolePhysics instantiates for " << mpl->GetParticleName() << " at " << mpl << " Mass "
94  << mass / CLHEP::GeV << " GeV Mag " << magn << " Process manager " << pmanager << G4endl;
95  }
96 
97  if (magn != 0.0) {
98  G4int idxt(0);
99  pmanager->RemoveProcess(idxt);
100  pmanager->AddProcess(new MonopoleTransportation(mpl, verbose), -1, 0, 0);
101  }
102 
103  if (mpl->GetPDGCharge() != 0.0) {
104  if (multiSc) {
105  G4hMultipleScattering* hmsc = new G4hMultipleScattering();
106  ph->RegisterProcess(hmsc, mpl);
107  }
108  G4hIonisation* hioni = new G4hIonisation();
109  ph->RegisterProcess(hioni, mpl);
110  }
111  if (magn != 0.0) {
112  G4mplIonisation* mplioni = new G4mplIonisation(magn);
113  if (!deltaRay) {
114  G4mplIonisationWithDeltaModel* ion = new G4mplIonisationWithDeltaModel(magn, "PAI");
115  ion->SetParticle(mpl);
116  mplioni->AddEmModel(0, ion, ion);
117  }
118  ph->RegisterProcess(mplioni, mpl);
119  }
120  pmanager->AddDiscreteProcess(new G4StepLimiter());
121  if (verbose > 1) {
122  pmanager->DumpInfo();
123  }
124  }
125  }
126 }
T getUntrackedParameter(std::string const &, T const &) const
~CMSMonopolePhysics() override
std::vector< int > pdgEncodings
HepPDT::ParticleDataTable ParticleDataTable
int ii
Definition: cuy.py:589
std::vector< Monopole * > monopoles
std::vector< int > elCharges
HepPDT::ParticleData ParticleData
std::vector< std::string > names
void ConstructProcess() override
void ConstructParticle() override
CMSMonopolePhysics(const HepPDT::ParticleDataTable *table_, const edm::ParameterSet &p)
G4double MagneticCharge() const
Definition: Monopole.h:20
std::vector< double > masses